Line data Source code
1 : module cam_history_support
2 :
3 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
4 : !!
5 : !! cam_history_support is used by cam_history as well as by the dycores
6 : !! (for vertical coordinate and "mdim" support). Some parameters are
7 : !! also referenced by cam_grid_support (although those could be copied
8 : !! if necessary).
9 : !!
10 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
11 :
12 : use shr_kind_mod, only: r8=>shr_kind_r8, shr_kind_cl, shr_kind_cxx
13 : use pio, only: var_desc_t, file_desc_t
14 : use cam_abortutils, only: endrun
15 : use cam_logfile, only: iulog
16 : use spmd_utils, only: masterproc
17 : use cam_grid_support, only: cam_grid_patch_t, cam_grid_header_info_t
18 : use cam_grid_support, only: max_hcoordname_len, maxsplitfiles
19 : use cam_pio_utils, only: cam_pio_handle_error
20 :
21 : implicit none
22 : private
23 : save
24 :
25 : integer, parameter, public :: max_string_len = shr_kind_cxx
26 : integer, parameter, public :: max_chars = shr_kind_cl ! max chars for char variables
27 : integer, parameter, public :: field_op_len = 3 ! max chars for field operation string (sum/dif)
28 : integer, parameter, public :: fieldname_len = 32 ! max chars for field name
29 : integer, parameter, public :: fieldname_suffix_len = 3 ! length of field name suffix ("&IC")
30 : integer, parameter, public :: fieldname_lenp2 = fieldname_len + 2 ! allow for extra characters
31 : ! max_fieldname_len = max chars for field name (including suffix)
32 : integer, parameter, public :: max_fieldname_len = fieldname_len + fieldname_suffix_len
33 :
34 : integer, parameter, public :: pflds = 1000 ! max number of fields for namelist entries fincl and fexcl
35 : ! also used in write restart
36 : integer, parameter, public :: ptapes = 12 ! max number of tapes
37 :
38 : ! A special symbol for declaring a field which has no vertical or
39 : ! non-grid dimensions. It is here (rather than cam_history) so that it
40 : ! be checked by add_hist_coord
41 : character(len=10), parameter, public :: horiz_only = 'horiz_only'
42 :
43 : type, public :: history_patch_t
44 : character(len=max_chars) :: namelist_entry = ''
45 : ! lon_axis_name and lat_axis_name are not used if collected_output = .true.
46 : character(len=max_fieldname_len) :: lon_axis_name = ''
47 : character(len=max_fieldname_len) :: lat_axis_name = ''
48 : logical :: collected_output
49 : ! There is one patch for every grid and one header_info for every patch
50 : type(cam_grid_patch_t), pointer :: patches(:) => NULL()
51 : type (cam_grid_header_info_t), pointer :: header_info(:) => NULL()
52 : contains
53 : procedure :: write_attrs => history_patch_write_attrs
54 : procedure :: write_vals => history_patch_write_vals
55 : procedure :: field_name => history_patch_field_name
56 : procedure :: num_hdims => history_patch_num_hdims
57 : procedure :: get_var_data => history_patch_get_var_data
58 : procedure :: write_var => history_patch_write_var
59 : procedure :: compact => history_patch_compact
60 : procedure :: active_cols => history_patch_active_cols
61 : procedure :: deallocate => history_patch_deallocate
62 : end type history_patch_t
63 :
64 : !
65 : ! dim_index_2d, dim_index_3d: 2-D & 3-D dimension index lower & upper bounds
66 : !
67 : type, public :: dim_index_2d ! 2-D dimension index
68 : integer :: beg1, end1 ! lower & upper bounds of 1st dimension
69 : integer :: beg2, end2 ! lower & upper bounds of 2nd dimension
70 : contains
71 : procedure :: dim_sizes_2d => dim_index_2d_dim_sizes_2d
72 : procedure :: dim_sizes_arr => dim_index_2d_dim_size_arr
73 : generic :: dim_sizes => dim_sizes_arr, dim_sizes_2d
74 : end type dim_index_2d
75 :
76 : type, public :: dim_index_3d ! 3-D dimension index
77 : integer :: beg1, end1 ! lower & upper bounds of 1st dimension
78 : integer :: beg2, end2 ! lower & upper bounds of 2nd dimension
79 : integer :: beg3, end3 ! lower & upper bounds of 3rd dimension
80 : contains
81 : procedure :: dim_sizes_3d => dim_index_3d_dim_sizes_3d
82 : procedure :: dim_sizes_arr => dim_index_3d_dim_size_arr
83 : generic :: dim_sizes => dim_sizes_arr, dim_sizes_3d
84 : end type dim_index_3d
85 :
86 : !---------------------------------------------------------------------------
87 : !
88 : ! field_info: A derived type containing information in an addfld call.
89 : !
90 : !---------------------------------------------------------------------------
91 : type, public :: field_info
92 :
93 : logical :: flag_xyfill ! non-applicable xy points flagged with fillvalue
94 : logical :: is_subcol ! .true. iff field output as subcol
95 : integer, pointer :: mdims(:) => NULL() ! indicies into hist_coords list
96 : integer, pointer :: shape(:) => NULL() ! shape of field on file
97 :
98 : real(r8) :: fillvalue ! fillvalue for this variable, set to default if not explicit in addfld
99 :
100 : integer :: numlev ! vertical dimension (.nc file and internal arr)
101 :
102 : integer :: begdim1 ! on-node dim1 start index
103 : integer :: enddim1 ! on-node dim1 end index
104 :
105 : integer :: begdim2 ! on-node dim2 start index
106 : integer :: enddim2 ! on-node dim2 end index
107 :
108 : integer :: begdim3 ! on-node chunk or lat start index
109 : integer :: enddim3 ! on-node chunk or lat end index
110 :
111 : logical :: colperchunk ! .true. iff ncols /= chunksize
112 :
113 : integer :: decomp_type ! type of decomposition (e.g., physics or dynamics)
114 :
115 : ! meridional and zonal complements are for fields defined as part of a
116 : ! 2-D vector. These IDs are used to facilitate interpolated history output
117 : ! At most one of these will be a positive field ID.
118 : integer :: meridional_complement ! meridional field id or -1
119 : integer :: zonal_complement ! zonal field id or -1
120 :
121 : ! Logical to determine if subcycle averages are allowed
122 : logical :: sampled_on_subcycle = .false.
123 :
124 : character(len=field_op_len) :: field_op = '' ! 'sum' or 'dif'
125 : integer :: op_field1_id ! first field id or -1
126 : integer :: op_field2_id ! second field id or -1
127 :
128 : character(len=max_fieldname_len) :: name ! field name
129 : character(len=max_chars) :: long_name ! long name
130 : character(len=max_chars) :: units ! units
131 : character(len=3) :: mixing_ratio ! 'dry' or 'wet'
132 : character(len=max_chars) :: sampling_seq ! sampling sequence - if not every timestep, how often field is sampled
133 : ! (i.e., how often "outfld" is called): every other; only during LW/SW
134 : ! radiation calcs; etc.
135 : character(len=max_chars) :: cell_methods ! optional cell_methods attribute
136 : contains
137 : procedure :: is_composed => field_info_is_composed
138 : procedure :: get_shape => field_info_get_shape
139 : procedure :: get_bounds => field_info_get_bounds
140 : procedure :: get_dims_2d => field_info_get_dims_2d
141 : procedure :: get_dims_3d => field_info_get_dims_3d
142 : generic :: get_dims => get_dims_2d, get_dims_3d
143 : end type field_info
144 :
145 : real(r8), parameter, public :: fillvalue = 1.e36_r8 ! fill value for netcdf fields
146 :
147 :
148 : !---------------------------------------------------------------------------
149 : !
150 : ! hentry: elements of an entry in the list of active fields on a single
151 : ! history file
152 : ! nacs is an accumulation counter which normally counts an entire
153 : ! chunk (physics) or block (dynamics) as accumulated as a single
154 : ! entity. The per-chunk counting avoids counting multiple outfld
155 : ! calls as multiple accumulations. Only the value of the first chunk
156 : ! or block is written to or read from a history restart file.
157 : ! For certain actions (e.g., only accumulating on
158 : ! non-fillvalue or accumulating based on local time), nacs has an
159 : ! entry for every column.
160 : ! nacs does not keep track of levels
161 : !
162 : !---------------------------------------------------------------------------
163 : type, public:: hentry
164 : type (field_info) :: field ! field information
165 : character(len=1) :: avgflag ! averaging flag
166 : character(len=max_chars) :: time_op ! time operator (e.g. max, min, avg)
167 : character(len=max_fieldname_len) :: op_field1 ! field1 name for sum or dif operation
168 : character(len=max_fieldname_len) :: op_field2 ! field2 name for sum or dif operation
169 :
170 : integer :: hwrt_prec ! history output precision
171 : real(r8), pointer :: hbuf(:,:,:) => NULL()
172 : real(r8), private :: hbuf_integral ! area weighted integral of active field
173 : real(r8), pointer :: sbuf(:,:,:) => NULL() ! for standard deviation
174 : real(r8), pointer :: wbuf(:,:) => NULL() ! pointer to area weights
175 : type(var_desc_t), pointer :: varid(:) => NULL() ! variable ids
176 : integer, pointer :: nacs(:,:) => NULL() ! accumulation counter
177 : type(var_desc_t), pointer :: nacs_varid => NULL()
178 : integer :: beg_nstep ! starting time step for nstep normalization
179 : type(var_desc_t), pointer :: beg_nstep_varid=> NULL()
180 : type(var_desc_t), pointer :: sbuf_varid => NULL()
181 : type(var_desc_t), pointer :: wbuf_varid => NULL()
182 : contains
183 : procedure :: get_global => hentry_get_global
184 : procedure :: put_global => hentry_put_global
185 : end type hentry
186 :
187 : !---------------------------------------------------------------------------
188 : !
189 : ! active_entry: derived type containing information for a history tape
190 : !
191 : !---------------------------------------------------------------------------
192 : type, public:: active_entry
193 :
194 : type(hentry), pointer :: hlist(:)
195 :
196 : integer, pointer :: grid_ids(:) => NULL()
197 : type(history_patch_t), pointer :: patches(:) => NULL()
198 :
199 : !
200 : ! PIO ids
201 : !
202 :
203 : type(file_desc_t) :: Files(maxsplitfiles) ! PIO file ids
204 :
205 : type(var_desc_t) :: mdtid ! var id for timestep
206 : type(var_desc_t) :: ndbaseid ! var id for base day
207 : type(var_desc_t) :: nsbaseid ! var id for base seconds of base day
208 : type(var_desc_t) :: nbdateid ! var id for base date
209 : type(var_desc_t) :: nbsecid ! var id for base seconds of base date
210 : type(var_desc_t) :: ndcurid ! var id for current day
211 : type(var_desc_t) :: nscurid ! var id for current seconds of current day
212 : type(var_desc_t) :: dateid ! var id for current date
213 : type(var_desc_t) :: co2vmrid ! var id for co2 volume mixing ratio
214 : type(var_desc_t) :: ch4vmrid ! var id for ch4 volume mixing ratio
215 : type(var_desc_t) :: n2ovmrid ! var id for n2o volume mixing ratio
216 : type(var_desc_t) :: f11vmrid ! var id for f11 volume mixing ratio
217 : type(var_desc_t) :: f12vmrid ! var id for f12 volume mixing ratio
218 : type(var_desc_t) :: sol_tsiid ! var id for total solar irradiance (W/m2)
219 : type(var_desc_t) :: datesecid ! var id for curent seconds of current date
220 : #if ( defined BFB_CAM_SCAM_IOP )
221 : type(var_desc_t) :: bdateid ! var id for base date
222 : type(var_desc_t) :: tsecid ! var id for curent seconds of current date
223 : #endif
224 : type(var_desc_t) :: nstephid ! var id for current timestep
225 : type(var_desc_t) :: timeid ! var id for time
226 : type(var_desc_t) :: tbndid ! var id for time_bounds
227 : type(var_desc_t) :: date_writtenid ! var id for date time sample written
228 : type(var_desc_t) :: time_writtenid ! var id for time time sample written
229 : type(var_desc_t) :: f107id ! var id for f107
230 : type(var_desc_t) :: f107aid ! var id for f107a
231 : type(var_desc_t) :: f107pid ! var id for f107p
232 : type(var_desc_t) :: kpid ! var id for kp
233 : type(var_desc_t) :: apid ! var id for ap
234 : type(var_desc_t) :: byimfid ! var id IMF BY
235 : type(var_desc_t) :: bzimfid ! var id IMF BZ
236 : type(var_desc_t) :: swvelid ! var id solar wind velocity
237 : type(var_desc_t) :: swdenid ! var id solar wind density
238 : type(var_desc_t) :: colat_crit1_id ! critical colatitude
239 : type(var_desc_t) :: colat_crit2_id ! critical colatitude
240 :
241 : end type active_entry
242 :
243 : !---------------------------------------------------------------------------
244 : !
245 : ! formula_terms_t: Information for formula terms (CF convention) variables
246 : ! Used to add a formula-terms variable to the history file
247 : ! Also adds a string, '<name>: <var_name>' to the parent
248 : ! mdim's 'formula_terms' attribute.
249 : !
250 : !---------------------------------------------------------------------------
251 : type, public :: formula_terms_t
252 : character(len=max_fieldname_len) :: a_name = '' ! 'A' term variable name
253 : character(len=max_string_len) :: a_long_name = '' ! 'A' long name
254 : real(r8), pointer :: a_values(:) => null() ! 'A' variable values
255 : character(len=max_fieldname_len) :: b_name = '' ! 'B' term variable name
256 : character(len=max_string_len) :: b_long_name = '' ! 'B' long name
257 : real(r8), pointer :: b_values(:) => null() ! 'B' variable values
258 : character(len=max_fieldname_len) :: p0_name = '' ! 'p0' term variable name
259 : character(len=max_string_len) :: p0_long_name = '' ! 'p0' long name
260 : character(len=max_chars) :: p0_units = '' ! 'p0' variable units
261 : real(r8) :: p0_value = fillvalue ! 'p0' variable values
262 : character(len=max_fieldname_len) :: ps_name = '' ! 'ps' term variable name
263 : end type formula_terms_t
264 :
265 : !---------------------------------------------------------------------------
266 : !
267 : ! hist_coord_t: Information for history variable dimension attributes
268 : !
269 : !---------------------------------------------------------------------------
270 : type, public :: hist_coord_t
271 : character(len=max_hcoordname_len) :: name = '' ! coordinate name
272 : integer :: dimsize = 0 ! size of dimension
273 : character(len=max_hcoordname_len) :: dimname = '' ! optional dimension name
274 : character(len=max_chars) :: long_name = '' ! 'long_name' attribute
275 : character(len=max_chars) :: units = '' ! 'units' attribute
276 : character(len=max_chars) :: bounds_name = '' ! 'bounds' attribute (& name of bounds variable)
277 : character(len=max_chars) :: standard_name = ''! 'standard_name' attribute
278 : character(len=4) :: positive = '' ! 'positive' attribute ('up' or 'down')
279 : integer, pointer :: integer_values(:) => null() ! dim values if integral
280 : real(r8), pointer :: real_values(:) => null() ! dim values if real
281 : real(r8), pointer :: bounds(:,:) => null() ! dim bounds
282 : type(formula_terms_t) :: formula_terms ! vars for formula terms
283 : logical :: integer_dim ! .true. iff dim has integral values
284 : logical :: vertical_coord ! .true. iff dim is vertical
285 : end type hist_coord_t
286 :
287 : ! Some parameters for use with interpolated output namelist items
288 : integer, parameter, public :: interp_type_native = 0
289 : integer, parameter, public :: interp_type_bilinear = 1
290 : integer, parameter, public :: interp_gridtype_equal_poles = 1
291 : integer, parameter, public :: interp_gridtype_gauss = 2
292 : integer, parameter, public :: interp_gridtype_equal_nopoles = 3
293 :
294 : !---------------------------------------------------------------------------
295 : !
296 : ! interp_info_t: Information for lat/lon interpolated history output
297 : !
298 : !---------------------------------------------------------------------------
299 : type, public :: interp_info_t
300 : ! store the lat-lon grid information
301 : character(len=28) :: gridname = ''
302 : integer :: grid_id = -1
303 : ! gridtype = 1 equally spaced, including poles (FV scalars output grid)
304 : ! gridtype = 2 Gauss grid (not implemented)
305 : ! gridtype = 3 equally spaced, no poles (FV staggered velocity)
306 : integer :: interp_gridtype = interp_gridtype_equal_poles
307 : ! interpolate_type = 0: native high order interpolation
308 : ! interpolate_type = 1: bilinear interpolation
309 : integer :: interp_type = interp_type_bilinear
310 : integer :: interp_nlat = 0
311 : integer :: interp_nlon = 0
312 : real(r8), pointer :: interp_lat(:) => NULL()
313 : real(r8), pointer :: interp_lon(:) => NULL()
314 : real(r8), pointer :: interp_gweight(:) => NULL()
315 : end type interp_info_t
316 :
317 : !! Coordinate variables
318 : integer, public :: registeredmdims = 0
319 : integer, public :: maxvarmdims = 1
320 : character(len=9), parameter, public :: mdim_var_name = 'mdimnames'
321 : integer, parameter :: maxmdims = 25 ! arbitrary limit
322 : type(hist_coord_t), public :: hist_coords(maxmdims)
323 :
324 : public :: add_hist_coord, add_vert_coord
325 : public :: write_hist_coord_attrs, write_hist_coord_vars
326 : public :: lookup_hist_coord_indices, hist_coord_find_levels
327 : public :: get_hist_coord_index, hist_coord_name, hist_coord_size
328 : public :: hist_dimension_name
329 : public :: hist_dimension_values
330 :
331 : interface add_hist_coord
332 : module procedure add_hist_coord_regonly
333 : module procedure add_hist_coord_int
334 : module procedure add_hist_coord_r8
335 : end interface
336 :
337 : interface hist_coord_size
338 : module procedure hist_coord_size_char
339 : module procedure hist_coord_size_int
340 : end interface hist_coord_size
341 :
342 : interface hist_dimension_values
343 : module procedure hist_dimension_values_r8
344 : module procedure hist_dimension_values_int
345 : end interface hist_dimension_values
346 :
347 : interface assignment(=)
348 : module procedure field_copy
349 : module procedure formula_terms_copy
350 : end interface
351 :
352 : interface check_hist_coord
353 : ! NB: This is supposed to be a private interface
354 : ! check_hist_coord: returns 0 if <name> is not registered as an mdim
355 : ! returns i if <name> is registered with compatible values
356 : ! calls endrun if <name> is registered with incompatible values
357 : ! Versions without the <name> argument return .true. or .false.
358 : module procedure check_hist_coord_char
359 : module procedure check_hist_coord_int
360 : module procedure check_hist_coord_int1
361 : module procedure check_hist_coord_r8
362 : module procedure check_hist_coord_r81
363 : module procedure check_hist_coord_r82
364 : module procedure check_hist_coord_ft
365 : module procedure check_hist_coord_all
366 : end interface
367 :
368 : !!---------------------------------------------------------------------------
369 :
370 : contains
371 :
372 76680960 : subroutine dim_index_2d_dim_sizes_2d(this, dim1, dim2)
373 :
374 : ! Dummy arguments
375 : class(dim_index_2d) :: this
376 : integer, intent(out) :: dim1
377 : integer, intent(out) :: dim2
378 :
379 76680960 : dim1 = MAX(0, this%end1 - this%beg1 + 1)
380 76680960 : dim2 = MAX(0, this%end2 - this%beg2 + 1)
381 :
382 76680960 : end subroutine dim_index_2d_dim_sizes_2d
383 :
384 0 : subroutine dim_index_2d_dim_size_arr(this, dims)
385 :
386 : ! Dummy arguments
387 : class(dim_index_2d) :: this
388 : integer, intent(out) :: dims(:)
389 :
390 0 : if (size(dims) < 2) then
391 0 : call endrun('dim_index_2d_dim_sizes: dims must have at least two elements')
392 : end if
393 :
394 0 : call this%dim_sizes(dims(1), dims(2))
395 :
396 0 : end subroutine dim_index_2d_dim_size_arr
397 :
398 4672512 : subroutine dim_index_3d_dim_sizes_3d(this, dim1, dim2, dim3)
399 :
400 : ! Dummy arguments
401 : class(dim_index_3d) :: this
402 : integer, intent(out) :: dim1
403 : integer, intent(out) :: dim2
404 : integer, intent(out) :: dim3
405 :
406 4672512 : dim1 = MAX(0, this%end1 - this%beg1 + 1)
407 4672512 : dim2 = MAX(0, this%end2 - this%beg2 + 1)
408 4672512 : dim3 = MAX(0, this%end3 - this%beg3 + 1)
409 :
410 4672512 : end subroutine dim_index_3d_dim_sizes_3d
411 :
412 4672512 : subroutine dim_index_3d_dim_size_arr(this, dims)
413 :
414 : ! Dummy arguments
415 : class(dim_index_3d) :: this
416 : integer, intent(out) :: dims(:)
417 :
418 4672512 : if (size(dims) < 3) then
419 0 : call endrun('dim_index_3d_dim_sizes: dims must have at least three elements')
420 : end if
421 :
422 4672512 : call this%dim_sizes(dims(1), dims(2), dims(3))
423 :
424 4672512 : end subroutine dim_index_3d_dim_size_arr
425 :
426 : ! field_info_get_dims_2d: Retrieve bounds for stepping through a chunk
427 155118336 : type(dim_index_2d) function field_info_get_dims_2d(this, col) result(dims)
428 : use cam_grid_support, only: cam_grid_get_block_count
429 :
430 : ! Dummy argument
431 : class(field_info) :: this
432 : integer, intent(in) :: col
433 :
434 : ! Local variable
435 : integer :: endi
436 :
437 155118336 : if (this%colperchunk) then
438 154556160 : endi = this%begdim1 + cam_grid_get_block_count(this%decomp_type, col) - 1
439 154556160 : dims = dim_index_2d(this%begdim1, endi, this%begdim2, this%enddim2)
440 : else
441 562176 : dims = dim_index_2d(this%begdim1, this%enddim1, this%begdim2, this%enddim2)
442 : end if
443 155118336 : end function field_info_get_dims_2d
444 :
445 : ! field_info_get_dims_3d: Retrieve grid decomp bounds
446 4672512 : type(dim_index_3d) function field_info_get_dims_3d(this) result(dims)
447 :
448 : ! Dummy argument
449 : class(field_info) :: this
450 :
451 : dims = dim_index_3d(this%begdim1, this%enddim1, this%begdim2, this%enddim2,&
452 4672512 : this%begdim3, this%enddim3)
453 :
454 155118336 : end function field_info_get_dims_3d
455 :
456 : ! field_info_is_composed: Return whether this field is composed of two other fields
457 9345024 : pure logical function field_info_is_composed(this)
458 : class(field_info), intent(IN) :: this
459 :
460 : field_info_is_composed = ((trim(adjustl(this%field_op))=='sum' .or. trim(adjustl(this%field_op))=='dif') .and. &
461 9345024 : this%op_field1_id /= -1 .and. this%op_field2_id /= -1)
462 9345024 : end function field_info_is_composed
463 :
464 : ! field_info_get_shape: Return a pointer to the field's global shape.
465 : ! Calculate it first if necessary
466 4672512 : subroutine field_info_get_shape(this, shape_out, rank_out)
467 : use cam_grid_support, only: cam_grid_dimensions
468 :
469 : ! Dummy arguments
470 : class(field_info) :: this
471 : integer, intent(out) :: shape_out(:)
472 : integer, intent(out) :: rank_out
473 :
474 : ! Local arguments
475 : integer :: rank, i, pos
476 : integer :: gdims(2)
477 :
478 4672512 : if (.not. associated(this%shape)) then
479 : ! Calculate field's global shape
480 1557504 : call cam_grid_dimensions(this%decomp_type, gdims, rank)
481 1557504 : pos = rank
482 1557504 : if (associated(this%mdims)) then
483 1195776 : rank = rank + size(this%mdims)
484 : end if
485 4672512 : allocate(this%shape(rank))
486 4672512 : this%shape(1:pos) = gdims(1:pos)
487 3115008 : if (rank > pos) then
488 1668096 : do i = 1, size(this%mdims)
489 834048 : pos = pos + 1
490 1668096 : this%shape(pos) = hist_coords(this%mdims(i))%dimsize
491 : end do
492 : end if
493 : end if
494 :
495 4672512 : rank_out = size(this%shape)
496 :
497 4672512 : if (size(shape_out) < rank_out) then
498 0 : call endrun('field_info_get_shape: shape_out too small')
499 : end if
500 :
501 16519680 : shape_out(1:rank_out) = this%shape(1:rank_out)
502 4672512 : if (size(shape_out) > rank_out) then
503 34877952 : shape_out(rank_out+1:) = 1
504 : end if
505 :
506 4672512 : end subroutine field_info_get_shape
507 :
508 15762432 : subroutine field_info_get_bounds(this, dim, beg, end)
509 :
510 : ! Dummy arguments
511 : class(field_info) :: this
512 : integer, intent(in) :: dim
513 : integer, intent(out) :: beg
514 : integer, intent(out) :: end
515 :
516 15762432 : select case(dim)
517 : case (1)
518 0 : beg = this%begdim1
519 0 : end = this%enddim1
520 : case (2)
521 0 : beg = this%begdim2
522 0 : end = this%enddim2
523 : case (3)
524 15762432 : beg = this%begdim3
525 15762432 : end = this%enddim3
526 : case default
527 15762432 : call endrun('field_info_get_bounds: dim must be 1, 2, or 3')
528 : end select
529 :
530 4672512 : end subroutine field_info_get_bounds
531 :
532 1932288 : subroutine hentry_get_global(this, gval)
533 :
534 : ! Dummy arguments
535 : class(hentry) :: this
536 : real(r8), intent(out) :: gval
537 :
538 1932288 : gval=this%hbuf_integral
539 :
540 1932288 : end subroutine hentry_get_global
541 :
542 966144 : subroutine hentry_put_global(this, gval)
543 :
544 : ! Dummy arguments
545 : class(hentry) :: this
546 : real(r8), intent(in) :: gval
547 :
548 966144 : this%hbuf_integral=gval
549 :
550 966144 : end subroutine hentry_put_global
551 :
552 : ! history_patch_write_attrs: Define coordinate variables and attributes
553 : ! for a patch
554 0 : subroutine history_patch_write_attrs(this, File)
555 : use cam_grid_support, only: cam_grid_is_unstructured
556 : use pio, only: file_desc_t, var_desc_t, pio_put_att, pio_double
557 : use cam_pio_utils, only: cam_pio_def_dim, cam_pio_def_var, cam_pio_handle_error
558 :
559 : ! Dummy arguments
560 : class(history_patch_t) :: this
561 : type(file_desc_t), intent(inout) :: File ! PIO file Handle
562 :
563 : ! Local variable
564 : type(cam_grid_patch_t), pointer :: patchptr
565 : type(var_desc_t), pointer :: vardesc_lat => NULL()
566 : type(var_desc_t), pointer :: vardesc_lon => NULL()
567 : character(len=128) :: errormsg
568 : character(len=max_chars) :: lat_name
569 : character(len=max_chars) :: lon_name
570 : character(len=max_chars) :: col_name
571 : character(len=max_chars) :: temp_str
572 : integer :: dimid1, dimid2 ! PIO dim IDs
573 : integer :: num_patches
574 : integer :: temp1, temp2
575 : integer :: latid, lonid ! Coordinate dims
576 : integer :: i, ierr
577 : logical :: col_only
578 : logical :: unstruct
579 : character(len=*), parameter :: subname = 'history_patch_write_attrs'
580 :
581 0 : num_patches = size(this%patches)
582 0 : if (associated(this%header_info)) then
583 : ! Make sure header_info is the right size
584 0 : if (size(this%header_info) /= num_patches) then
585 0 : write(errormsg, '(a,2(i0,a))') 'Size mismatch between header_info (', &
586 0 : size(this%header_info), ') and patches (', num_patches, ')'
587 0 : call endrun(subname//': '//errormsg)
588 : end if
589 : else
590 0 : allocate(this%header_info(num_patches))
591 : end if
592 :
593 : ! Write attributes for each patch
594 0 : do i = 1, num_patches
595 0 : patchptr => this%patches(i)
596 0 : call this%header_info(i)%set_gridid(patchptr%gridid())
597 0 : unstruct = cam_grid_is_unstructured(patchptr%gridid())
598 : ! What are the dimension(s) for this patch?
599 0 : col_only = this%collected_output
600 0 : if (num_patches == 1) then
601 : ! Backwards compatibility
602 0 : if (unstruct .or. col_only) then
603 0 : col_name = 'ncol'
604 : else
605 0 : col_name = ''
606 : end if
607 0 : lat_name = 'lat'
608 0 : lon_name = 'lon'
609 : else
610 0 : call patchptr%get_axis_names(lat_name, lon_name, col_name, col_only)
611 : end if
612 : ! Define the dimensions (latx/lonx or ncolx)
613 : ! col_name is set for unstructured output or collected columns (ncolx)
614 0 : if (len_trim(col_name) > 0) then
615 0 : call patchptr%get_global_size(gsize=temp1)
616 0 : if (temp1 <= 0) then
617 0 : call endrun(subname//': col dimsize must be positive')
618 : end if
619 0 : if (unstruct .and. (.not. col_only)) then
620 : ! For the case of unstructured output without collected column
621 : ! output, we need to make sure that the ncolx dimension is unique
622 0 : col_name = trim(col_name)//'_'//trim(this%lon_axis_name)//'_'//trim(this%lat_axis_name)
623 : end if
624 0 : call cam_pio_def_dim(File, trim(col_name), temp1, dimid1, existOK=.false.)
625 0 : call this%header_info(i)%set_hdims(dimid1)
626 0 : latid = dimid1
627 0 : lonid = dimid1
628 : else
629 0 : lat_name = trim(lat_name)//'_'//trim(this%lat_axis_name)
630 0 : call patchptr%get_global_size(temp1, temp2)
631 0 : if (temp1 <= 0) then
632 0 : call endrun(subname//': lat dimsize must be positive')
633 : end if
634 0 : call cam_pio_def_dim(File, trim(lat_name), temp1, dimid1, existOK=.true.)
635 0 : latid = dimid1
636 0 : lon_name = trim(lon_name)//'_'//trim(this%lon_axis_name)
637 0 : if (temp2 <= 0) then
638 0 : call endrun(subname//': lon dimsize must be positive')
639 : end if
640 0 : call cam_pio_def_dim(File, trim(lon_name), temp2, dimid2, existOK=.true.)
641 0 : lonid = dimid2
642 0 : call this%header_info(i)%set_hdims(lonid, latid)
643 : end if
644 : !! Define the latx (coordinate) variable
645 0 : if (unstruct .and. (.not. col_only)) then
646 : ! We need to make sure the latx name is unique
647 0 : lat_name = trim(lat_name)//'_'//trim(this%lon_axis_name)//'_'//trim(this%lat_axis_name)
648 : end if
649 0 : allocate(vardesc_lat)
650 : call cam_pio_def_var(File, trim(lat_name), pio_double, (/latid/), &
651 0 : vardesc_lat, existOK=.true.)
652 : ! Coordinate attributes
653 0 : call patchptr%get_coord_long_name('lat', temp_str)
654 0 : if (len_trim(temp_str) > 0) then
655 0 : ierr = pio_put_att(File, vardesc_lat, 'long_name', trim(temp_str))
656 0 : call cam_pio_handle_error(ierr, subname//': Unable to define long_name')
657 : end if
658 0 : call patchptr%get_coord_units('lat', temp_str)
659 0 : if (len_trim(temp_str) > 0) then
660 0 : ierr = pio_put_att(File, vardesc_lat, 'units', trim(temp_str))
661 0 : call cam_pio_handle_error(ierr, subname//': Unable to define units')
662 : end if
663 : !! Define the lonx (coordinate) variable
664 0 : if (unstruct .and. (.not. col_only)) then
665 : ! We need to make sure the lonx name is unique
666 0 : lon_name = trim(lon_name)//'_'//trim(this%lon_axis_name)//'_'//trim(this%lat_axis_name)
667 : end if
668 0 : allocate(vardesc_lon)
669 : call cam_pio_def_var(File, trim(lon_name), pio_double, (/lonid/), &
670 0 : vardesc_lon, existOK=.true.)
671 : ! Coordinate attributes
672 0 : call patchptr%get_coord_long_name('lon', temp_str)
673 0 : if (len_trim(temp_str) > 0) then
674 0 : ierr = pio_put_att(File, vardesc_lon, 'long_name', trim(temp_str))
675 0 : call cam_pio_handle_error(ierr, subname//': Unable to define long_name')
676 : end if
677 0 : call patchptr%get_coord_units('lon', temp_str)
678 0 : if (len_trim(temp_str) > 0) then
679 0 : ierr = pio_put_att(File, vardesc_lon, 'units', trim(temp_str))
680 0 : call cam_pio_handle_error(ierr, subname//': Unable to define units')
681 : end if
682 0 : call this%header_info(i)%set_varids(vardesc_lon, vardesc_lat)
683 0 : nullify(vardesc_lat, vardesc_lon) ! They belong to the header_info now
684 : end do
685 :
686 0 : end subroutine history_patch_write_attrs
687 :
688 : ! history_patch_write_vals: Write coordinate variable values for a patch
689 0 : subroutine history_patch_write_vals(this, File)
690 0 : use pio, only: file_desc_t, var_desc_t
691 :
692 : ! Dummy arguments
693 : class(history_patch_t) :: this
694 : type(file_desc_t), intent(inout) :: File ! PIO file Handle
695 :
696 : ! Local variable
697 : type(cam_grid_patch_t), pointer :: patchptr
698 : type(var_desc_t), pointer :: vardesc => NULL() ! PIO var desc
699 : character(len=128) :: errormsg
700 : integer :: num_patches
701 : integer :: i
702 :
703 0 : num_patches = size(this%patches)
704 0 : if (.not. associated(this%header_info)) then
705 : ! We need this for dim and coord variable descriptors
706 0 : write(errormsg, '(2a)') 'No header info for ', trim(this%namelist_entry)
707 0 : call endrun('history_patch_write_vals: '//errormsg)
708 : end if
709 :
710 : ! Write attributes for each patch
711 0 : do i = 1, num_patches
712 0 : patchptr => this%patches(i)
713 : ! Write the coordinate variables (or just lat/lon for column output)
714 0 : call patchptr%write_coord_vals(File, this%header_info(i))
715 : end do
716 :
717 0 : end subroutine history_patch_write_vals
718 :
719 : ! history_patch_field_name: Add patch description to field name
720 0 : subroutine history_patch_field_name(this, name)
721 : ! Dummy arguments
722 : class(history_patch_t) :: this
723 : character(len=*), intent(inout) :: name
724 :
725 0 : if (.not. this%collected_output) then
726 : ! Add patch description info to the variable name
727 0 : name = trim(name)//'_'//trim(this%lon_axis_name)//'_'//trim(this%lat_axis_name)
728 : end if
729 0 : end subroutine history_patch_field_name
730 :
731 : ! history_patch_num_hdims: Find the number of horizontal dimensions for
732 : ! the indicated grid
733 0 : integer function history_patch_num_hdims(this, gridid)
734 : ! Dummy arguments
735 : class(history_patch_t) :: this
736 : integer, intent(in) :: gridid ! The field's grid
737 :
738 : ! Local variables
739 : type(cam_grid_patch_t), pointer :: patchptr
740 : character(len=128) :: errormsg
741 : integer :: i
742 : integer :: num_patches
743 :
744 : ! Basic sanity checks, is this patch OK?
745 0 : num_patches = size(this%patches)
746 0 : if (associated(this%header_info)) then
747 : ! Make sure header_info is the right size
748 0 : if (size(this%header_info) /= num_patches) then
749 0 : write(errormsg, '(a,2(i0,a))') 'Size mismatch between header_info (', &
750 0 : size(this%header_info), ') and patches (', num_patches, ')'
751 0 : call endrun('history_patch_num_hdims: '//errormsg)
752 : end if
753 : else
754 0 : write(errormsg, *) 'No header info for patch, ', trim(this%namelist_entry)
755 0 : call endrun('history_patch_num_hdims: '//errormsg)
756 : end if
757 :
758 : ! Find the correct patch by matching grid ID
759 0 : history_patch_num_hdims = -1
760 0 : do i = 1, num_patches
761 0 : patchptr => this%patches(i)
762 0 : if (patchptr%gridid() == gridid) then
763 : ! This is the right patch, set the return val and quit loop
764 0 : history_patch_num_hdims = this%header_info(i)%num_hdims()
765 0 : exit
766 0 : else if (i >= num_patches) then
767 0 : write(errormsg, '(3a,i0)') 'No grid found for patch, ', &
768 0 : trim(this%namelist_entry), '. Was looking for decomp ', gridid
769 0 : call endrun('history_patch_num_hdims: '//errormsg)
770 : ! No else needed
771 : end if
772 : end do
773 0 : if (history_patch_num_hdims <= 0) then
774 0 : write(errormsg, '(2a,2(a,i0))') 'INTERNAL: No grid patch for ', &
775 0 : trim(this%namelist_entry), ', num_patches = ',num_patches, &
776 0 : ', gridid = ', gridid
777 0 : call endrun('history_patch_num_hdims: '//errormsg)
778 : end if
779 :
780 0 : end function history_patch_num_hdims
781 :
782 : ! history_patch_get_var_data: Calculate data relevant to history variable
783 : ! on a patch by substituting patch dimension ids for the horiz. ids
784 : ! and adding patch information to the variable name
785 0 : subroutine history_patch_get_var_data(this, name, dimids, gridid)
786 : ! Dummy arguments
787 : class(history_patch_t) :: this
788 : character(len=*), intent(inout) :: name
789 : integer, intent(inout) :: dimids(:) ! Grid dimids
790 : integer, intent(in) :: gridid ! The field's grid
791 :
792 : ! Local variables
793 : type(cam_grid_patch_t), pointer :: patchptr
794 : type (cam_grid_header_info_t), pointer :: histptr
795 : character(len=128) :: errormsg
796 : integer :: num_patches
797 : integer :: i
798 :
799 : ! Basic sanity checks, is this patch OK?
800 0 : num_patches = size(this%patches)
801 0 : if (associated(this%header_info)) then
802 : ! Make sure header_info is the right size
803 0 : if (size(this%header_info) /= num_patches) then
804 0 : write(errormsg, '(a,2(i0,a))') 'Size mismatch between header_info (', &
805 0 : size(this%header_info), ') and patches (', num_patches, ')'
806 0 : call endrun('history_patch_get_var_data: '//errormsg)
807 : end if
808 : else
809 0 : write(errormsg, *) 'No header info for patch, ', trim(this%namelist_entry)
810 0 : call endrun('history_patch_get_var_data: '//errormsg)
811 : end if
812 :
813 : ! Find the correct patch by matching grid ID
814 0 : do i = 1, num_patches
815 0 : patchptr => this%patches(i)
816 0 : if (patchptr%gridid() == gridid) then
817 : ! This is the right patch, quit loop
818 0 : histptr => this%header_info(i)
819 0 : exit
820 0 : else if (i >= num_patches) then
821 0 : write(errormsg, '(3a,i0)') 'No grid found for patch, ', &
822 0 : trim(this%namelist_entry), '. Was looking for decomp ', gridid
823 0 : call endrun('history_patch_get_var_data: '//errormsg)
824 : ! No else needed
825 : end if
826 : end do
827 :
828 : ! We have the correct patch, replace the horizontal dimension ids
829 0 : do i = 1, histptr%num_hdims()
830 0 : dimids(i) = histptr%get_hdimid(i)
831 : end do
832 : ! Re-define the variable name
833 0 : call this%field_name(name)
834 :
835 0 : end subroutine history_patch_get_var_data
836 :
837 0 : subroutine history_patch_compact(this)
838 :
839 : ! Dummy arguments
840 : class(history_patch_t) :: this
841 :
842 : ! Local variables
843 : integer :: num_patches
844 : integer :: i
845 :
846 0 : num_patches = size(this%patches)
847 :
848 : ! Find the correct patch by matching grid ID
849 0 : do i = 1, num_patches
850 0 : call this%patches(i)%compact(this%collected_output)
851 : end do
852 :
853 0 : end subroutine history_patch_compact
854 :
855 0 : subroutine history_patch_write_var(this, File, gridid, adims, dtype, hbuf, varid)
856 : use pio, only: file_desc_t, var_desc_t, io_desc_t
857 : use pio, only: pio_write_darray
858 : use cam_pio_utils, only: cam_pio_handle_error, cam_pio_var_info
859 :
860 : ! Dummy arguments
861 : class(history_patch_t) :: this
862 : type(file_desc_t), intent(inout) :: File ! PIO file handle
863 : integer, intent(in) :: gridid
864 : integer, intent(in) :: adims(:)
865 : integer, intent(in) :: dtype
866 : real(r8), intent(in) :: hbuf(:,:,:)
867 : type(var_desc_t), pointer :: varid
868 :
869 : ! Local variables
870 : type(cam_grid_patch_t), pointer :: patchptr
871 : character(len=128) :: errormsg
872 : integer :: num_patches
873 : integer :: i, idx
874 : integer :: uid ! unlimited dim ID
875 : type(io_desc_t), pointer :: iodesc
876 : integer :: ierr, nfdims
877 : integer :: fdimlens(7), dimids(7)
878 :
879 0 : num_patches = size(this%patches)
880 :
881 : ! Find the correct patch by matching grid ID
882 0 : do i = 1, num_patches
883 0 : patchptr => this%patches(i)
884 0 : if (patchptr%gridid() == gridid) then
885 : ! This is the right patch, quit loop
886 : exit
887 0 : else if (i >= num_patches) then
888 0 : write(errormsg, '(3a,i0)') 'No grid found for patch, ', &
889 0 : trim(this%namelist_entry), '. Was looking for decomp ', gridid
890 0 : call endrun('history_patch_write_var: '//trim(errormsg))
891 : ! No else needed
892 : end if
893 : end do
894 :
895 : ! We have the right grid, write the hbuf
896 0 : call cam_pio_var_info(File, varid, nfdims, dimids, fdimlens, unlimDimID=uid)
897 0 : idx = 1
898 0 : do i = 1, nfdims
899 0 : if (i > idx) then
900 0 : dimids(idx) = dimids(i)
901 : end if
902 0 : if (dimids(i) /= uid) then
903 0 : idx = idx + 1
904 : end if
905 : end do
906 0 : nfdims = nfdims - COUNT(dimids(1:nfdims) == uid)
907 0 : call patchptr%get_decomp(adims, fdimlens(1:nfdims), dtype, iodesc)
908 0 : if (size(adims) == 2) then
909 0 : call pio_write_darray(File, varid, iodesc, hbuf(:,1,:), ierr)
910 0 : else if (size(adims) == 3) then
911 0 : call pio_write_darray(File, varid, iodesc, hbuf, ierr)
912 : else
913 0 : call endrun("history_patch_write_var: adims must be rank 2 or 3")
914 : end if
915 0 : call cam_pio_handle_error(ierr, 'history_patch_write_var: Error writing variable')
916 :
917 0 : end subroutine history_patch_write_var
918 :
919 0 : subroutine history_patch_active_cols(this, gridid, lchnk, active)
920 : ! Dummy arguments
921 : class(history_patch_t) :: this
922 : integer, intent(in) :: gridid ! desired grid
923 : integer, intent(in) :: lchnk ! chunk or block number
924 : logical, intent(out) :: active(:)
925 :
926 : ! Local variables
927 : type(cam_grid_patch_t), pointer :: patchptr
928 : character(len=128) :: errormsg
929 : integer :: num_patches
930 : integer :: i
931 :
932 0 : num_patches = size(this%patches)
933 :
934 : ! Find the correct patch by matching grid ID
935 0 : do i = 1, num_patches
936 0 : patchptr => this%patches(i)
937 0 : if (patchptr%gridid() == gridid) then
938 : ! This is the right patch, quit loop
939 : exit
940 0 : else if (i >= num_patches) then
941 0 : write(errormsg, '(3a,i0)') 'No grid found for patch, ', &
942 0 : trim(this%namelist_entry), '. Was looking for decomp ', gridid
943 0 : call endrun('history_patch_active_cols: '//errormsg)
944 : ! No else needed
945 : end if
946 : end do
947 :
948 : ! If we get here, patchptr is the grid patch we want
949 0 : call patchptr%active_cols(lchnk, active)
950 :
951 0 : end subroutine history_patch_active_cols
952 :
953 0 : subroutine history_patch_deallocate(this)
954 : ! Dummy argument
955 : class(history_patch_t) :: this
956 : ! Local variable
957 : integer :: i
958 :
959 0 : this%lon_axis_name = ''
960 0 : this%lat_axis_name = ''
961 :
962 0 : if (associated(this%patches)) then
963 0 : do i = 1, size(this%patches)
964 0 : call this%patches(i)%deallocate()
965 : end do
966 0 : deallocate(this%patches)
967 0 : nullify(this%patches)
968 : end if
969 :
970 0 : if (associated(this%header_info)) then
971 0 : do i = 1, size(this%header_info)
972 0 : call this%header_info(i)%deallocate()
973 : end do
974 0 : deallocate(this%header_info)
975 0 : nullify(this%header_info)
976 : end if
977 :
978 0 : end subroutine history_patch_deallocate
979 :
980 : subroutine field_copy(f_out, f_in)
981 : type(field_info), intent(in) :: f_in
982 : type(field_info), intent(out) :: f_out
983 :
984 : f_out%flag_xyfill= f_in%flag_xyfill
985 : f_out%is_subcol = f_in%is_subcol
986 : f_out%fillvalue= f_in%fillvalue
987 : f_out%numlev = f_in%numlev ! vertical dimension (.nc file and internal arr)
988 : f_out%begdim1 = f_in%begdim1 ! on-node dim1 start index
989 : f_out%enddim1 = f_in%enddim1 ! on-node dim1 end index
990 : f_out%begdim2 = f_in%begdim2 ! on-node dim2 start index
991 : f_out%enddim2 = f_in%enddim2 ! on-node dim2 end index
992 : f_out%begdim3 = f_in%begdim3 ! on-node chunk or lat start index
993 : f_out%enddim3 = f_in%enddim3 ! on-node chunk or lat end index
994 : f_out%decomp_type = f_in%decomp_type ! type of decomposition (physics or dynamics)
995 :
996 : f_out%meridional_complement = f_in%meridional_complement ! id or -1
997 : f_out%zonal_complement = f_in%zonal_complement ! id or -1
998 : f_out%field_op = f_in%field_op ! sum,dif, or ''
999 : f_out%op_field1_id = f_in%op_field1_id ! id or -1
1000 : f_out%op_field2_id = f_in%op_field2_id ! id or -1
1001 :
1002 : f_out%name = f_in%name ! field name
1003 : f_out%long_name = f_in%long_name ! long name
1004 : f_out%units = f_in%units ! units
1005 : f_out%mixing_ratio = f_in%mixing_ratio ! mixing_ratio
1006 : f_out%sampling_seq = f_in%sampling_seq ! sampling sequence - if not every timestep, how often field is sampled
1007 : f_out%cell_methods = f_in%cell_methods
1008 :
1009 : if(associated(f_in%mdims)) then
1010 : f_out%mdims=>f_in%mdims
1011 : else
1012 : nullify(f_out%mdims)
1013 : end if
1014 :
1015 : end subroutine field_copy
1016 :
1017 3072 : subroutine formula_terms_copy(f_out, f_in)
1018 : type(formula_terms_t), intent(in) :: f_in
1019 : type(formula_terms_t), intent(out) :: f_out
1020 :
1021 3072 : f_out%a_name = f_in%a_name
1022 3072 : f_out%a_long_name = f_in%a_long_name
1023 3072 : f_out%a_values => f_in%a_values
1024 3072 : f_out%b_name = f_in%b_name
1025 3072 : f_out%b_long_name = f_in%b_long_name
1026 3072 : f_out%b_values => f_in%b_values
1027 3072 : f_out%p0_name = f_in%p0_name
1028 3072 : f_out%p0_long_name = f_in%p0_long_name
1029 3072 : f_out%p0_units = f_in%p0_units
1030 3072 : f_out%p0_value = f_in%p0_value
1031 3072 : f_out%ps_name = f_in%ps_name
1032 3072 : end subroutine formula_terms_copy
1033 :
1034 9518592 : integer function get_hist_coord_index(mdimname)
1035 : ! Input variables
1036 : character(len=*), intent(in) :: mdimname
1037 : ! Local variable
1038 : integer :: i
1039 :
1040 9518592 : get_hist_coord_index = -1
1041 9681408 : do i = 1, registeredmdims
1042 9681408 : if(trim(mdimname) == trim(hist_coords(i)%name)) then
1043 9500160 : get_hist_coord_index = i
1044 9500160 : exit
1045 : end if
1046 : end do
1047 :
1048 9518592 : end function get_hist_coord_index
1049 :
1050 6144 : character(len=max_hcoordname_len) function hist_coord_name(index)
1051 : ! Input variables
1052 : integer, intent(in) :: index
1053 :
1054 6144 : if ((index > 0) .and. (index <= registeredmdims)) then
1055 6144 : hist_coord_name = hist_coords(index)%name
1056 : else
1057 0 : call endrun('hist_coord_name: index out of range')
1058 : end if
1059 :
1060 6144 : end function hist_coord_name
1061 :
1062 9490944 : integer function hist_coord_size_int(index)
1063 : ! Input variables
1064 : integer, intent(in) :: index
1065 :
1066 9490944 : if (index > 0) then
1067 9490944 : hist_coord_size_int = hist_coords(index)%dimsize
1068 : else
1069 : hist_coord_size_int = -1
1070 : end if
1071 :
1072 9490944 : end function hist_coord_size_int
1073 :
1074 0 : integer function hist_coord_size_char(mdimname)
1075 : ! Input variables
1076 : character(len=*), intent(in) :: mdimname
1077 : ! Local variable
1078 : integer :: i
1079 :
1080 0 : i = get_hist_coord_index(mdimname)
1081 0 : hist_coord_size_char = hist_coord_size(i)
1082 :
1083 0 : end function hist_coord_size_char
1084 :
1085 : ! Functions to check consistent term definition for hist coords
1086 0 : logical function check_hist_coord_char(defined, input)
1087 :
1088 : ! Input variables
1089 : character(len=*), intent(in) :: defined
1090 : character(len=*), intent(in), optional :: input
1091 :
1092 0 : if (len_trim(defined) == 0) then
1093 : ! In this case, we assume the current value is undefined so any input OK
1094 : check_hist_coord_char = .true.
1095 0 : else if (present(input)) then
1096 : ! We have to match definitions
1097 0 : check_hist_coord_char = (trim(input) == trim(defined))
1098 : else
1099 : ! Not sure here. We have a value and are redefining without one?
1100 : check_hist_coord_char = .false.
1101 : end if
1102 0 : end function check_hist_coord_char
1103 :
1104 0 : logical function check_hist_coord_int(defined, input)
1105 :
1106 : ! Input variables
1107 : integer, intent(in) :: defined
1108 : integer, intent(in), optional :: input
1109 :
1110 0 : if (defined == 0) then
1111 : ! In this case, we assume the current value is undefined so any input OK
1112 : check_hist_coord_int = .true.
1113 0 : else if (present(input)) then
1114 : ! We have to match definitions
1115 0 : check_hist_coord_int = (input == defined)
1116 : else
1117 : ! Not sure here. We have a value and are redefining without one?
1118 : check_hist_coord_int = .false.
1119 : end if
1120 0 : end function check_hist_coord_int
1121 :
1122 0 : logical function check_hist_coord_int1(defined, input)
1123 :
1124 : ! Input variables
1125 : integer, pointer :: defined(:)
1126 : integer, intent(in), optional :: input(:)
1127 :
1128 : ! Local variables
1129 : integer :: i
1130 :
1131 0 : if (.not. associated(defined)) then
1132 : ! In this case, we assume the current value is undefined so any input OK
1133 : check_hist_coord_int1 = .true.
1134 0 : else if (present(input)) then
1135 : ! We have to match definitions
1136 0 : check_hist_coord_int1 = (size(input) == size(defined))
1137 : else
1138 : ! Not sure here. We have a value and are redefining without one?
1139 : check_hist_coord_int1 = .false.
1140 : end if
1141 0 : if (check_hist_coord_int1 .and. associated(defined)) then
1142 : ! Need to check the values
1143 0 : do i = 1, size(defined)
1144 0 : if (defined(i) /= input(i)) then
1145 : check_hist_coord_int1 = .false.
1146 : exit
1147 : end if
1148 : end do
1149 : end if
1150 0 : end function check_hist_coord_int1
1151 :
1152 0 : logical function check_hist_coord_r8(defined, input)
1153 :
1154 : ! Input variables
1155 : real(r8), intent(in) :: defined
1156 : real(r8), intent(in), optional :: input
1157 :
1158 0 : if (defined == fillvalue) then
1159 : ! In this case, we assume the current value is undefined so any input OK
1160 : check_hist_coord_r8 = .true.
1161 0 : else if (present(input)) then
1162 : ! We have to match definitions
1163 0 : check_hist_coord_r8 = (input == defined)
1164 : else
1165 : ! Not sure here. We have a value and are redefining without one?
1166 : check_hist_coord_r8 = .false.
1167 : end if
1168 0 : end function check_hist_coord_r8
1169 :
1170 0 : logical function check_hist_coord_r81(defined, input)
1171 :
1172 : ! Input variables
1173 : real(r8), pointer :: defined(:)
1174 : real(r8), intent(in), optional :: input(:)
1175 :
1176 : ! Local variables
1177 : integer :: i
1178 :
1179 0 : if (.not. associated(defined)) then
1180 : ! In this case, we assume the current value is undefined so any input OK
1181 : check_hist_coord_r81 = .true.
1182 0 : else if (present(input)) then
1183 : ! We have to match definitions
1184 0 : check_hist_coord_r81 = (size(input) == size(defined))
1185 : else
1186 : ! Not sure here. We have a value and are redefining without one?
1187 : check_hist_coord_r81 = .false.
1188 : end if
1189 0 : if (check_hist_coord_r81 .and. associated(defined)) then
1190 : ! Need to check the values
1191 0 : do i = 1, size(defined)
1192 0 : if (defined(i) /= input(i)) then
1193 : check_hist_coord_r81 = .false.
1194 : exit
1195 : end if
1196 : end do
1197 : end if
1198 0 : end function check_hist_coord_r81
1199 :
1200 0 : logical function check_hist_coord_r82(defined, input)
1201 :
1202 : ! Input variables
1203 : real(r8), pointer :: defined(:,:)
1204 : real(r8), intent(in), optional :: input(:,:)
1205 :
1206 : ! Local variables
1207 : integer :: i, j
1208 :
1209 0 : if (.not. associated(defined)) then
1210 : ! In this case, we assume the current value is undefined so any input OK
1211 : check_hist_coord_r82 = .true.
1212 0 : else if (present(input)) then
1213 : ! We have to match definitions
1214 : check_hist_coord_r82 = ((size(input, 1) == size(defined, 1)) .and. &
1215 0 : (size(input, 2) == size(defined, 2)))
1216 : else
1217 : ! Not sure here. We have a value and are redefining without one?
1218 : check_hist_coord_r82 = .false.
1219 : end if
1220 0 : if (check_hist_coord_r82 .and. associated(defined)) then
1221 : ! Need to check the values
1222 0 : do j = 1, size(defined, 2)
1223 0 : do i = 1, size(defined, 1)
1224 0 : if (defined(i, j) /= input(i, j)) then
1225 : check_hist_coord_r82 = .false.
1226 : exit
1227 : end if
1228 : end do
1229 : end do
1230 : end if
1231 0 : end function check_hist_coord_r82
1232 :
1233 0 : logical function check_hist_coord_ft(defined, input)
1234 :
1235 : ! Input variables
1236 : type(formula_terms_t), intent(in) :: defined
1237 : type(formula_terms_t), intent(in), optional :: input
1238 :
1239 : ! We will assume that if formula_terms has been defined, a_name has a value
1240 0 : if (len_trim(defined%a_name) == 0) then
1241 : ! In this case, we assume the current value is undefined so any input OK
1242 : check_hist_coord_ft = .true.
1243 0 : else if (present(input)) then
1244 : ! We have to match definitions
1245 : ! Need to check the values
1246 : check_hist_coord_ft = &
1247 : check_hist_coord(defined%a_name, input%a_name) .and. &
1248 : check_hist_coord(defined%a_long_name, input%a_long_name) .and. &
1249 : check_hist_coord(defined%a_values, input%a_values) .and. &
1250 : check_hist_coord(defined%b_name, input%b_name) .and. &
1251 : check_hist_coord(defined%b_long_name, input%b_long_name) .and. &
1252 : check_hist_coord(defined%b_values, input%b_values) .and. &
1253 : check_hist_coord(defined%p0_name, input%p0_name) .and. &
1254 : check_hist_coord(defined%p0_long_name, input%p0_long_name) .and. &
1255 : check_hist_coord(defined%p0_units, input%p0_units) .and. &
1256 : check_hist_coord(defined%p0_value, input%p0_value) .and. &
1257 0 : check_hist_coord(defined%ps_name, input%ps_name)
1258 : else
1259 : ! Not sure here. We have a value and are redefining without one?
1260 : check_hist_coord_ft = .false.
1261 : end if
1262 0 : end function check_hist_coord_ft
1263 :
1264 : ! check_hist_coord: returns 0 if <name> is not registered as a hist coord
1265 : ! returns i if <name> is registered with compatible values
1266 : ! calls endrun if <name> is registered with incompatible values
1267 12288 : integer function check_hist_coord_all(name, vlen, long_name, units, bounds, &
1268 12288 : i_values, r_values, bounds_name, positive, standard_name, formula_terms)
1269 :
1270 : ! Input variables
1271 : character(len=*), intent(in) :: name
1272 : integer, intent(in) :: vlen
1273 : character(len=*), intent(in), optional :: long_name
1274 : character(len=*), intent(in), optional :: units
1275 : character(len=*), intent(in), optional :: bounds_name
1276 : integer, intent(in), optional :: i_values(:)
1277 : real(r8), intent(in), optional :: r_values(:)
1278 : real(r8), intent(in), optional :: bounds(:,:)
1279 : character(len=*), intent(in), optional :: positive
1280 : character(len=*), intent(in), optional :: standard_name
1281 : type(formula_terms_t), intent(in), optional :: formula_terms
1282 :
1283 : ! Local variables
1284 : character(len=120) :: errormsg
1285 : integer :: i
1286 :
1287 12288 : i = get_hist_coord_index(trim(name))
1288 : ! If i > 0, this mdim has already been registered
1289 12288 : if (i > 0) then
1290 0 : check_hist_coord_all = i
1291 0 : if (.not. check_hist_coord(hist_coords(i)%dimsize, vlen)) then
1292 0 : write(errormsg, *) 'ERROR: Attempt to register dimension, '//trim(name)//' with incompatible size'
1293 0 : call endrun(errormsg)
1294 : end if
1295 0 : if (.not. check_hist_coord(hist_coords(i)%long_name, long_name)) then
1296 0 : write(errormsg, *) 'ERROR: Attempt to register dimension, ',trim(name),' with a different long_name'
1297 0 : call endrun(errormsg)
1298 : end if
1299 0 : if (.not. check_hist_coord(hist_coords(i)%units, units)) then
1300 0 : write(errormsg, *) 'ERROR: Attempt to register dimension, ',trim(name),' with different units'
1301 0 : call endrun(errormsg)
1302 : end if
1303 0 : if (.not. check_hist_coord(hist_coords(i)%bounds_name, bounds_name)) then
1304 0 : write(errormsg, *) 'ERROR: Attempt to register dimension, ',trim(name),' with a different bounds_name'
1305 0 : call endrun(errormsg)
1306 : end if
1307 0 : if (.not. check_hist_coord(hist_coords(i)%standard_name, standard_name)) then
1308 0 : write(errormsg, *) 'ERROR: Attempt to register dimension, ',trim(name),' with a different standard_name'
1309 0 : call endrun(errormsg)
1310 : end if
1311 0 : if (.not. check_hist_coord(hist_coords(i)%positive, positive)) then
1312 0 : write(errormsg, *) 'ERROR: Attempt to register dimension, ',trim(name),' with a different value of positive'
1313 0 : call endrun(errormsg)
1314 : end if
1315 : ! Since the integer_dim defaults to .true., double check which to check
1316 0 : if ((.not. hist_coords(i)%integer_dim) .or. &
1317 : associated(hist_coords(i)%real_values)) then
1318 0 : if (.not. check_hist_coord(hist_coords(i)%real_values, r_values)) then
1319 0 : write(errormsg, *) 'ERROR: Attempt to register dimension, ',trim(name),' with different values'
1320 0 : call endrun(errormsg)
1321 0 : else if (present(i_values)) then
1322 0 : write(errormsg, *) 'ERROR: Attempt to register integer values for real dimension'
1323 0 : call endrun(errormsg)
1324 : end if
1325 : else
1326 0 : if (.not. check_hist_coord(hist_coords(i)%integer_values, i_values)) then
1327 0 : write(errormsg, *) 'ERROR: Attempt to register dimension, ',trim(name),' with different values'
1328 0 : call endrun(errormsg)
1329 0 : else if (present(i_values) .and. present(r_values)) then
1330 0 : write(errormsg, *) 'ERROR: Attempt to register real values for integer dimension'
1331 0 : call endrun(errormsg)
1332 : end if
1333 : end if
1334 0 : if (.not. check_hist_coord(hist_coords(i)%bounds, bounds)) then
1335 0 : write(errormsg, *) 'ERROR: Attempt to register dimension, ',trim(name),' with different bounds'
1336 0 : call endrun(errormsg)
1337 : end if
1338 0 : if (.not. check_hist_coord(hist_coords(i)%formula_terms, formula_terms)) then
1339 0 : write(errormsg, *) 'ERROR: Attempt to register dimension, ',trim(name),' with different formula_terms'
1340 0 : call endrun(errormsg)
1341 : end if
1342 : else
1343 : check_hist_coord_all = 0
1344 : end if
1345 12288 : end function check_hist_coord_all
1346 :
1347 6144 : subroutine add_hist_coord_regonly(name, index)
1348 :
1349 : ! Input variable
1350 : character(len=*), intent(in) :: name
1351 : integer, optional, intent(out) :: index
1352 :
1353 : ! Local variables
1354 : character(len=120) :: errormsg
1355 : integer :: i
1356 :
1357 6144 : if ((trim(name) == trim(horiz_only)) .or. (len_trim(name) == 0)) then
1358 0 : call endrun('ADD_HIST_COORD: '//trim(name)//' is not a valid coordinate name')
1359 : end if
1360 6144 : i = get_hist_coord_index(trim(name))
1361 : ! If i > 0, this mdim has already been registered
1362 6144 : if (i <= 0) then
1363 6144 : registeredmdims = registeredmdims + 1
1364 6144 : if (registeredmdims > maxmdims) then
1365 0 : call endrun('Too many dimensions in add_hist_coord.')
1366 : end if
1367 6144 : if (len_trim(name) > max_hcoordname_len) then
1368 0 : write(errormsg,'(a,i3,a)') 'History coord name exceeds the ', &
1369 0 : max_hcoordname_len, ' character length limit'
1370 0 : call endrun(errormsg)
1371 : end if
1372 6144 : hist_coords(registeredmdims)%name = trim(name)
1373 6144 : hist_coords(registeredmdims)%dimsize = 0
1374 6144 : hist_coords(registeredmdims)%long_name = ''
1375 6144 : hist_coords(registeredmdims)%units = ''
1376 6144 : hist_coords(registeredmdims)%integer_dim = .true.
1377 6144 : hist_coords(registeredmdims)%positive = ''
1378 6144 : hist_coords(registeredmdims)%standard_name = ''
1379 6144 : if (present(index)) then
1380 6144 : index = registeredmdims
1381 : end if
1382 : else
1383 0 : if (present(index)) then
1384 0 : index = i
1385 : end if
1386 : end if
1387 :
1388 6144 : end subroutine add_hist_coord_regonly
1389 :
1390 0 : subroutine add_hist_coord_int(name, vlen, long_name, units, values, &
1391 : positive, standard_name, dimname)
1392 :
1393 : ! Input variables
1394 : character(len=*), intent(in) :: name
1395 : integer, intent(in) :: vlen
1396 : character(len=*), intent(in) :: long_name
1397 : character(len=*), intent(in), optional :: units
1398 : integer, intent(in), target, optional :: values(:)
1399 : character(len=*), intent(in), optional :: positive
1400 : character(len=*), intent(in), optional :: standard_name
1401 : character(len=*), intent(in), optional :: dimname
1402 :
1403 : ! Local variables
1404 : integer :: i
1405 :
1406 : ! First, check to see if it is OK to add this coord
1407 : i = check_hist_coord(name, vlen=vlen, long_name=long_name, units=units, &
1408 0 : i_values=values, positive=positive, standard_name=standard_name)
1409 : ! Register the name if necessary
1410 0 : if (i == 0) then
1411 0 : call add_hist_coord(trim(name), i)
1412 0 : if(masterproc) then
1413 0 : write(iulog, '(3a,i0,a,i0)') 'Registering hist coord: ', trim(name), &
1414 0 : '(', i, ') with length: ', vlen
1415 : end if
1416 : end if
1417 :
1418 : ! Set the coord's values
1419 0 : hist_coords(i)%dimsize = vlen
1420 0 : if (len_trim(long_name) > max_chars) then
1421 0 : if(masterproc) then
1422 0 : write(iulog,*) 'WARNING: long_name for ',trim(name),' too long'
1423 : end if
1424 : end if
1425 0 : hist_coords(i)%long_name = trim(long_name)
1426 0 : if (present(units)) then
1427 0 : hist_coords(i)%units = trim(units)
1428 : else
1429 0 : hist_coords(i)%units = ''
1430 : end if
1431 0 : hist_coords(i)%integer_dim = .true.
1432 0 : if (present(values)) then
1433 0 : hist_coords(i)%integer_values => values
1434 : endif
1435 0 : if (present(positive)) then
1436 0 : hist_coords(i)%positive = trim(positive)
1437 : end if
1438 0 : if (present(standard_name)) then
1439 0 : hist_coords(i)%standard_name = trim(standard_name)
1440 : end if
1441 0 : hist_coords(i)%vertical_coord = .false.
1442 0 : if (present(dimname)) then
1443 0 : hist_coords(i)%dimname = trim(dimname)
1444 : else
1445 0 : hist_coords(i)%dimname = ''
1446 : end if
1447 :
1448 0 : end subroutine add_hist_coord_int
1449 :
1450 6144 : subroutine add_hist_coord_r8(name, vlen, long_name, units, values, &
1451 6144 : bounds_name, bounds, positive, standard_name, vertical_coord, dimname)
1452 :
1453 : ! Input variables
1454 : character(len=*), intent(in) :: name
1455 : integer, intent(in) :: vlen
1456 : character(len=*), intent(in) :: long_name
1457 : character(len=*), intent(in) :: units
1458 : real(r8), intent(in), target :: values(:)
1459 : character(len=*), intent(in), optional :: bounds_name
1460 : real(r8), intent(in), target, optional :: bounds(:,:)
1461 : character(len=*), intent(in), optional :: positive
1462 : character(len=*), intent(in), optional :: standard_name
1463 : logical, intent(in), optional :: vertical_coord
1464 : character(len=*), intent(in), optional :: dimname
1465 :
1466 : ! Local variables
1467 : character(len=120) :: errormsg
1468 : integer :: i
1469 :
1470 : ! First, check to see if it is OK to add this coord
1471 : i = check_hist_coord(name, vlen=vlen, long_name=long_name, units=units, &
1472 : r_values=values, positive=positive, standard_name=standard_name, &
1473 21504 : bounds_name=bounds_name, bounds=bounds)
1474 : ! Register the name if necessary
1475 6144 : if (i == 0) then
1476 6144 : call add_hist_coord(trim(name), i)
1477 6144 : if(masterproc) then
1478 8 : write(iulog, '(3a,i0,a,i0)') 'Registering hist coord: ', trim(name), &
1479 16 : '(', i, ') with length: ', vlen
1480 : end if
1481 : end if
1482 :
1483 : ! Set the coord's size
1484 6144 : hist_coords(i)%dimsize = vlen
1485 6144 : if (len_trim(long_name) > max_chars) then
1486 0 : if(masterproc) then
1487 0 : write(iulog,*) 'WARNING: long_name for ',trim(name),' too long'
1488 : end if
1489 : end if
1490 6144 : hist_coords(i)%long_name = trim(long_name)
1491 6144 : if (len_trim(units) > 0) then
1492 6144 : hist_coords(i)%units = trim(units)
1493 : else
1494 0 : hist_coords(i)%units = '1'
1495 : end if
1496 6144 : hist_coords(i)%integer_dim = .false.
1497 6144 : hist_coords(i)%real_values => values
1498 6144 : if (present(positive)) then
1499 6144 : hist_coords(i)%positive = trim(positive)
1500 : end if
1501 6144 : if (present(standard_name)) then
1502 3072 : hist_coords(i)%standard_name = trim(standard_name)
1503 : end if
1504 6144 : if (present(bounds_name)) then
1505 0 : hist_coords(i)%bounds_name = trim(bounds_name)
1506 0 : if (.not. present(bounds)) then
1507 0 : write(errormsg,*) 'bounds must be present for ',trim(bounds_name)
1508 0 : call endrun(errormsg)
1509 : end if
1510 0 : hist_coords(i)%bounds => bounds
1511 6144 : else if (present(bounds)) then
1512 0 : write(errormsg,*) 'bounds_name must be present for bounds values'
1513 0 : call endrun(errormsg)
1514 : else
1515 6144 : hist_coords(i)%bounds_name = ''
1516 : end if
1517 6144 : if (present(vertical_coord)) then
1518 6144 : hist_coords(i)%vertical_coord = vertical_coord
1519 : else
1520 0 : hist_coords(i)%vertical_coord = .false.
1521 : end if
1522 6144 : if (present(dimname)) then
1523 0 : hist_coords(i)%dimname = trim(dimname)
1524 : else
1525 6144 : hist_coords(i)%dimname = ''
1526 : end if
1527 :
1528 6144 : end subroutine add_hist_coord_r8
1529 :
1530 6144 : subroutine add_vert_coord(name, vlen, long_name, units, values, &
1531 : positive, standard_name, formula_terms)
1532 :
1533 : ! Input variables
1534 : character(len=*), intent(in) :: name
1535 : integer, intent(in) :: vlen
1536 : character(len=*), intent(in) :: long_name
1537 : character(len=*), intent(in) :: units
1538 : real(r8), intent(in), target :: values(:)
1539 : character(len=*), intent(in), optional :: positive
1540 : character(len=*), intent(in), optional :: standard_name
1541 : type(formula_terms_t), intent(in), optional :: formula_terms
1542 :
1543 : ! Local variable
1544 : integer :: i
1545 :
1546 : ! First, check to see if it is OK to add this coord
1547 : i = check_hist_coord(name, vlen=vlen, long_name=long_name, units=units, &
1548 : r_values=values, positive=positive, standard_name=standard_name, &
1549 9216 : formula_terms=formula_terms)
1550 : ! Register the name and hist_coord values if necessary
1551 6144 : if (i == 0) then
1552 : call add_hist_coord(trim(name), vlen, long_name, units, values, &
1553 : positive=positive, standard_name=standard_name, &
1554 6144 : vertical_coord=.true.)
1555 6144 : i = get_hist_coord_index(trim(name))
1556 6144 : if(masterproc) then
1557 8 : write(iulog, '(3a,i0,a,i0)') 'Registering hist coord: ', trim(name), &
1558 16 : '(', i, ') with length: ', vlen
1559 : end if
1560 : end if
1561 :
1562 6144 : if (present(formula_terms)) then
1563 3072 : hist_coords(i)%formula_terms = formula_terms
1564 : end if
1565 :
1566 6144 : end subroutine add_vert_coord
1567 :
1568 36864 : subroutine write_hist_coord_attr(File, mdimind, boundsdim, dimonly, mdimid)
1569 : use pio, only: file_desc_t, var_desc_t, pio_put_att, pio_noerr, &
1570 : pio_int, pio_double, pio_inq_varid, pio_def_var
1571 : use cam_pio_utils, only: cam_pio_def_dim, cam_pio_def_var
1572 :
1573 : ! Input variables
1574 : type(file_desc_t), intent(inout) :: File ! PIO file Handle
1575 : integer, intent(in) :: mdimind ! Internal dim index
1576 : integer, intent(in) :: boundsdim ! Bounds dimension ID
1577 : logical, intent(in) :: dimonly ! No def_var if .true.
1578 : integer, optional, intent(out) :: mdimid
1579 :
1580 : ! Local variables
1581 : integer :: dimid ! PIO dimension ID
1582 : type(var_desc_t) :: vardesc ! PIO variable descriptor
1583 : character(len=120) :: errormsg
1584 : character(len=max_chars) :: formula_terms ! Constructed string
1585 : integer :: ierr
1586 : integer :: dtype
1587 : logical :: defvar ! True if var exists
1588 :
1589 : ! Create or check dimension for this coordinate
1590 36864 : if (len_trim(hist_coords(mdimind)%dimname) > 0) then
1591 : ! Dim can already exist if different from coord name
1592 : call cam_pio_def_dim(File, trim(hist_coords(mdimind)%dimname), &
1593 : hist_coords(mdimind)%dimsize, dimid, &
1594 : existOK=(trim(hist_coords(mdimind)%dimname) /= &
1595 0 : trim(hist_coords(mdimind)%name)))
1596 : else
1597 : ! The dimension has the same name as the coord -- must be new dim
1598 : call cam_pio_def_dim(File, trim(hist_coords(mdimind)%name), &
1599 36864 : hist_coords(mdimind)%dimsize, dimid, existOK=.false.)
1600 : end if
1601 : ! If the caller wants to know the NetCDF dimension ID, set it here
1602 36864 : if (present(mdimid)) then
1603 36864 : mdimid = dimid
1604 : end if
1605 36864 : if (.not. dimonly) then
1606 : ! Time to define the variable (only if there are values)
1607 36864 : if (hist_coords(mdimind)%integer_dim) then
1608 0 : dtype = pio_int
1609 0 : defvar = associated(hist_coords(mdimind)%integer_values)
1610 : else
1611 36864 : dtype = pio_double
1612 36864 : defvar = associated(hist_coords(mdimind)%real_values)
1613 : end if
1614 36864 : if (defvar) then
1615 : call cam_pio_def_var(File, trim(hist_coords(mdimind)%name), dtype, &
1616 73728 : (/dimid/), vardesc, existOK=.false.)
1617 : ! long_name
1618 36864 : ierr=pio_put_att(File, vardesc, 'long_name', trim(hist_coords(mdimind)%long_name))
1619 36864 : call cam_pio_handle_error(ierr, 'Error writing "long_name" attr in write_hist_coord_attr')
1620 : ! units
1621 36864 : if(len_trim(hist_coords(mdimind)%units) > 0) then
1622 : ierr=pio_put_att(File, vardesc, 'units', &
1623 36864 : trim(hist_coords(mdimind)%units))
1624 36864 : call cam_pio_handle_error(ierr, 'Error writing "units" attr in write_hist_coord_attr')
1625 : end if
1626 : ! positive
1627 36864 : if(len_trim(hist_coords(mdimind)%positive) > 0) then
1628 : ierr=pio_put_att(File, vardesc, 'positive', &
1629 36864 : trim(hist_coords(mdimind)%positive))
1630 36864 : call cam_pio_handle_error(ierr, 'Error writing "positive" attr in write_hist_coord_attr')
1631 : end if
1632 : ! standard_name
1633 36864 : if(len_trim(hist_coords(mdimind)%standard_name) > 0) then
1634 : ierr=pio_put_att(File, vardesc, 'standard_name', &
1635 18432 : trim(hist_coords(mdimind)%standard_name))
1636 18432 : call cam_pio_handle_error(ierr, 'Error writing "standard_name" attr in write_hist_coord_attr')
1637 : end if
1638 : ! formula_terms
1639 36864 : if(len_trim(hist_coords(mdimind)%formula_terms%a_name) > 0) then
1640 : write(formula_terms, '("a: ",a," b: ",a," p0: ",a," ps: ",a)') &
1641 18432 : trim(hist_coords(mdimind)%formula_terms%a_name), &
1642 18432 : trim(hist_coords(mdimind)%formula_terms%b_name), &
1643 18432 : trim(hist_coords(mdimind)%formula_terms%p0_name), &
1644 36864 : trim(hist_coords(mdimind)%formula_terms%ps_name)
1645 18432 : ierr=pio_put_att(File, vardesc, 'formula_terms', trim(formula_terms))
1646 18432 : call cam_pio_handle_error(ierr, 'Error writing "formula_terms" attr in write_hist_coord_attr')
1647 : end if
1648 : ! bounds
1649 36864 : if (associated(hist_coords(mdimind)%bounds)) then
1650 : ! Write name of the bounds variable
1651 0 : ierr=pio_put_att(File, vardesc, 'bounds', trim(hist_coords(mdimind)%bounds_name))
1652 0 : call cam_pio_handle_error(ierr, 'Error writing "bounds" attr in write_hist_coord_attr')
1653 : end if
1654 : end if
1655 :
1656 : ! Now, we need to define and populate the associated bounds variable
1657 : ! NB: Reusing vardesc, no longer assocated with main variable
1658 36864 : if (associated(hist_coords(mdimind)%bounds)) then
1659 0 : if (size(hist_coords(mdimind)%bounds,2) /= hist_coords(mdimind)%dimsize) then
1660 : ! If anyone hits this check, add a new dimension for this case
1661 0 : write(errormsg, *) 'The bounds variable, ', &
1662 0 : trim(hist_coords(mdimind)%bounds_name), &
1663 0 : ', needs to have dimension (2,', hist_coords(mdimind)%dimsize
1664 0 : call endrun(errormsg)
1665 : end if
1666 : call cam_pio_def_var(File, trim(hist_coords(mdimind)%bounds_name), &
1667 0 : pio_double, (/boundsdim,dimid/), vardesc, existOK=.false.)
1668 : end if
1669 :
1670 : ! See if we have formula_terms variables to define
1671 : ! Define the "a" variable name
1672 : ! NB: Reusing vardesc, no longer assocated with previous variables
1673 36864 : if (associated(hist_coords(mdimind)%formula_terms%a_values)) then
1674 18432 : if (size(hist_coords(mdimind)%formula_terms%a_values) /= hist_coords(mdimind)%dimsize) then
1675 0 : write(errormsg, *) 'The forumla_terms variable, ', &
1676 0 : trim(hist_coords(mdimind)%formula_terms%a_name), &
1677 0 : ', needs to have dimension', hist_coords(mdimind)%dimsize
1678 0 : call endrun(errormsg)
1679 : end if
1680 : call cam_pio_def_var(File, trim(hist_coords(mdimind)%formula_terms%a_name), &
1681 36864 : pio_double, (/dimid/), vardesc, existOK=.false.)
1682 18432 : ierr = pio_put_att(File, vardesc, 'long_name', trim(hist_coords(mdimind)%formula_terms%a_long_name))
1683 18432 : call cam_pio_handle_error(ierr, 'Error writing "long_name" attr for "a" formula_term in write_hist_coord_attr')
1684 : end if
1685 : ! Define the "b" variable name
1686 : ! NB: Reusing vardesc, no longer assocated with previous variables
1687 36864 : if (associated(hist_coords(mdimind)%formula_terms%b_values)) then
1688 18432 : if (size(hist_coords(mdimind)%formula_terms%b_values) /= hist_coords(mdimind)%dimsize) then
1689 0 : write(errormsg, *) 'The forumla_terms variable, ', &
1690 0 : trim(hist_coords(mdimind)%formula_terms%b_name), &
1691 0 : ', needs to have dimension', hist_coords(mdimind)%dimsize
1692 0 : call endrun(errormsg)
1693 : end if
1694 : call cam_pio_def_var(File, trim(hist_coords(mdimind)%formula_terms%b_name), &
1695 36864 : pio_double, (/dimid/), vardesc, existOK=.false.)
1696 18432 : ierr = pio_put_att(File, vardesc, 'long_name', trim(hist_coords(mdimind)%formula_terms%b_long_name))
1697 18432 : call cam_pio_handle_error(ierr, 'Error writing "long_name" attr for "b" formula_term in write_hist_coord_attr')
1698 : end if
1699 : ! Maybe define the p0 variable (this may be defined already which is OK)
1700 : ! NB: Reusing vardesc, no longer assocated with previous variables
1701 36864 : if (hist_coords(mdimind)%formula_terms%p0_value /= fillvalue) then
1702 18432 : ierr = pio_inq_varid(File, trim(hist_coords(mdimind)%formula_terms%p0_name), vardesc)
1703 18432 : if (ierr /= PIO_NOERR) then
1704 : ierr = pio_def_var(File, trim(hist_coords(mdimind)%formula_terms%p0_name), &
1705 9216 : pio_double, vardesc)
1706 9216 : call cam_pio_handle_error(ierr, 'Unable to define "p0" formula_terms variable in write_hist_coord_attr')
1707 9216 : ierr = pio_put_att(File, vardesc, 'long_name', trim(hist_coords(mdimind)%formula_terms%p0_long_name))
1708 9216 : call cam_pio_handle_error(ierr, 'Error writing "long_name" attr for "p0" formula_term in write_hist_coord_attr')
1709 9216 : ierr = pio_put_att(File, vardesc, 'units', trim(hist_coords(mdimind)%formula_terms%p0_units))
1710 9216 : call cam_pio_handle_error(ierr, 'Error writing "units" attr for "p0" formula_term in write_hist_coord_attr')
1711 : end if
1712 : end if
1713 : ! PS is not our responsibility
1714 : end if ! (.not. dimonly)
1715 :
1716 36864 : end subroutine write_hist_coord_attr
1717 :
1718 : !---------------------------------------------------------------------------
1719 : !
1720 : ! write_hist_coord_attrs
1721 : !
1722 : ! Write the dimension and coordinate attributes for the defined history
1723 : ! coordinates.
1724 : !
1725 : !---------------------------------------------------------------------------
1726 9216 : subroutine write_hist_coord_attrs(File, boundsdim, mdimids, writemdims_in)
1727 36864 : use pio, only: file_desc_t, var_desc_t, pio_put_att, &
1728 : pio_bcast_error, pio_internal_error, pio_seterrorhandling, &
1729 : pio_char
1730 : use cam_pio_utils, only: cam_pio_def_dim, cam_pio_def_var
1731 :
1732 : ! Input variables
1733 : type(file_desc_t), intent(inout) :: File ! PIO file Handle
1734 : integer, intent(in) :: boundsdim ! Bounds dimension ID
1735 : integer, optional, allocatable, intent(out) :: mdimids(:) ! NetCDF dim IDs
1736 : logical, optional, intent(in) :: writemdims_in ! Write mdim variable
1737 :
1738 : ! Local variables
1739 : integer :: i
1740 : integer :: ierr
1741 : integer :: dimids(2) ! PIO dimension IDs
1742 : logical :: writemdims ! Define an mdim variable
1743 : type(var_desc_t) :: vardesc ! PIO variable descriptor
1744 :
1745 9216 : if (present(mdimids)) then
1746 27648 : allocate(mdimids(registeredmdims))
1747 : end if
1748 :
1749 : ! We will handle errors for this routine
1750 9216 : call pio_seterrorhandling(File, PIO_BCAST_ERROR)
1751 :
1752 9216 : if (present(writemdims_in)) then
1753 9216 : writemdims = writemdims_in
1754 : else
1755 0 : writemdims = .false.
1756 : end if
1757 :
1758 : ! NB: Currently, writemdims is for restart and we don't need to write
1759 : ! these out in a history-restart file. This could change in the future.
1760 : ! which would require a change to the function of the fourth argument
1761 : ! Fill in the attribute information for each mdim
1762 46080 : do i = 1, registeredmdims
1763 46080 : if (present(mdimids)) then
1764 36864 : call write_hist_coord_attr(File, i, boundsdim, writemdims, mdimids(i))
1765 : else
1766 0 : call write_hist_coord_attr(File, i, boundsdim, writemdims)
1767 : end if
1768 : end do
1769 :
1770 9216 : if (writemdims) then
1771 : call cam_pio_def_dim(File, 'mdimslen', max_hcoordname_len, dimids(1), &
1772 0 : existOK=.true.)
1773 : call cam_pio_def_dim(File, 'num_mdims', registeredmdims, dimids(2), &
1774 0 : existOK=.true.)
1775 : call cam_pio_def_var(File, mdim_var_name, pio_char, dimids, vardesc, &
1776 0 : existOK=.false.)
1777 0 : ierr = pio_put_att(File, vardesc, 'long_name', 'mdim dimension names')
1778 0 : call cam_pio_handle_error(ierr, 'Error writing "long_name" attr for mdimnames in write_hist_coord_attrs')
1779 :
1780 : end if
1781 :
1782 : ! Back to I/O or die trying
1783 9216 : call pio_seterrorhandling(File, PIO_INTERNAL_ERROR)
1784 9216 : end subroutine write_hist_coord_attrs
1785 :
1786 36864 : subroutine write_hist_coord_var(File, mdimind)
1787 9216 : use pio, only: file_desc_t, var_desc_t, pio_put_var, pio_inq_varid
1788 :
1789 : ! Input variables
1790 : type(file_desc_t), intent(inout) :: File ! PIO file Handle
1791 : integer, intent(in) :: mdimind ! Internal dim index
1792 :
1793 : ! Local variables
1794 : type(var_desc_t) :: vardesc ! PIO variable descriptor
1795 : integer :: ierr
1796 :
1797 36864 : if ((hist_coords(mdimind)%integer_dim .and. &
1798 73728 : associated(hist_coords(mdimind)%integer_values)) .or. &
1799 : ((.not. hist_coords(mdimind)%integer_dim) .and. &
1800 : associated(hist_coords(mdimind)%real_values))) then
1801 : ! Check to make sure the variable already exists in the file
1802 36864 : ierr = pio_inq_varid(File, trim(hist_coords(mdimind)%name), vardesc)
1803 36864 : call cam_pio_handle_error(ierr, 'Error writing values for nonexistent dimension variable write_hist_coord_var')
1804 : ! Write out the values for this dimension variable
1805 36864 : if (hist_coords(mdimind)%integer_dim) then
1806 0 : ierr = pio_put_var(File, vardesc, hist_coords(mdimind)%integer_values)
1807 : else
1808 36864 : ierr = pio_put_var(File, vardesc, hist_coords(mdimind)%real_values)
1809 : end if
1810 36864 : call cam_pio_handle_error(ierr, 'Error writing variable values in write_hist_coord_var')
1811 : end if
1812 :
1813 : ! Now, we need to possibly write values for the associated bounds variable
1814 36864 : if (associated(hist_coords(mdimind)%bounds)) then
1815 : ! Check to make sure the variable already exists in the file
1816 : ! NB: Reusing vardesc, no longer assocated with previous variables
1817 0 : ierr = pio_inq_varid(File, trim(hist_coords(mdimind)%bounds_name), vardesc)
1818 0 : call cam_pio_handle_error(ierr, 'Error writing values for nonexistent bounds variable write_hist_coord_var')
1819 : ! Write out the values for this bounds variable
1820 0 : ierr = pio_put_var(File, vardesc, hist_coords(mdimind)%bounds)
1821 0 : call cam_pio_handle_error(ierr, 'Error writing bounds values in write_hist_coord_var')
1822 : end if
1823 :
1824 : ! Write values for the "a" variable name
1825 36864 : if (associated(hist_coords(mdimind)%formula_terms%a_values)) then
1826 : ! Check to make sure the variable already exists in the file
1827 : ! NB: Reusing vardesc, no longer assocated with previous variables
1828 18432 : ierr = pio_inq_varid(File, trim(hist_coords(mdimind)%formula_terms%a_name), vardesc)
1829 18432 : call cam_pio_handle_error(ierr, 'Error writing values for nonexistent "a" formula_terms variable write_hist_coord_var')
1830 : ! Write out the values for this "a" formula_terms variable
1831 18432 : ierr = pio_put_var(File, vardesc, hist_coords(mdimind)%formula_terms%a_values)
1832 18432 : call cam_pio_handle_error(ierr, 'Error writing "a" formula_terms values in write_hist_coord_var')
1833 : end if
1834 : ! Write values for the "b" variable name
1835 36864 : if (associated(hist_coords(mdimind)%formula_terms%b_values)) then
1836 : ! Check to make sure the variable already exists in the file
1837 : ! NB: Reusing vardesc, no longer assocated with previous variables
1838 18432 : ierr = pio_inq_varid(File, trim(hist_coords(mdimind)%formula_terms%b_name), vardesc)
1839 18432 : call cam_pio_handle_error(ierr, 'Error writing values for nonexistent "b" formula_terms variable write_hist_coord_var')
1840 : ! Write out the values for this "b" formula_terms variable
1841 18432 : ierr = pio_put_var(File, vardesc, hist_coords(mdimind)%formula_terms%b_values)
1842 18432 : call cam_pio_handle_error(ierr, 'Error writing "b" formula_terms values in write_hist_coord_var')
1843 : end if
1844 : ! Write values for the "p0" variable name (this may be an overwrite, too bad)
1845 36864 : if (hist_coords(mdimind)%formula_terms%p0_value /= fillvalue) then
1846 : ! Check to make sure the variable already exists in the file
1847 : ! NB: Reusing vardesc, no longer assocated with previous variables
1848 18432 : ierr = pio_inq_varid(File, trim(hist_coords(mdimind)%formula_terms%p0_name), vardesc)
1849 18432 : call cam_pio_handle_error(ierr, 'Error writing values for nonexistent "p0" formula_terms variable write_hist_coord_var')
1850 : ! Write out the values for this "p0" formula_terms variable
1851 18432 : ierr = pio_put_var(File, vardesc, hist_coords(mdimind)%formula_terms%p0_value)
1852 18432 : call cam_pio_handle_error(ierr, 'Error writing "p0" formula_terms values in write_hist_coord_var')
1853 : end if
1854 :
1855 36864 : end subroutine write_hist_coord_var
1856 :
1857 9216 : subroutine write_hist_coord_vars(File, writemdims_in)
1858 : use pio, only: file_desc_t, var_desc_t, pio_put_var, &
1859 : pio_bcast_error, pio_internal_error, &
1860 : pio_seterrorhandling, pio_inq_varid
1861 :
1862 : ! Input variables
1863 : type(file_desc_t), intent(inout) :: File ! PIO file Handle
1864 : logical, optional, intent(in) :: writemdims_in ! Write mdim variable
1865 :
1866 : ! Local variables
1867 : integer :: i
1868 : integer :: ierr
1869 : logical :: writemdims ! Define an mdim variable
1870 : type(var_desc_t) :: vardesc ! PIO variable descriptor
1871 9216 : character(len=max_hcoordname_len), allocatable :: mdimnames(:)
1872 :
1873 : ! We will handle errors for this routine
1874 9216 : call pio_seterrorhandling(File, PIO_BCAST_ERROR)
1875 :
1876 9216 : if (present(writemdims_in)) then
1877 9216 : writemdims = writemdims_in
1878 : else
1879 : writemdims = .false.
1880 : end if
1881 :
1882 9216 : if (writemdims) then
1883 0 : allocate(mdimnames(registeredmdims))
1884 : end if
1885 :
1886 : ! Write out the variable values for each mdim
1887 46080 : do i = 1, registeredmdims
1888 36864 : if (.not. writemdims) then
1889 : ! NB: Currently, writemdims is for restart and we don't need to write
1890 : ! these out in a history-restart file. This could change in the future
1891 : ! which is why it is a separate if block
1892 : ! Fill in the attribute information for each mdim
1893 36864 : call write_hist_coord_var(File, i)
1894 : end if
1895 9216 : if (writemdims) then
1896 0 : mdimnames(i) = trim(hist_coords(i)%name)
1897 : end if
1898 : end do
1899 :
1900 9216 : if (writemdims) then
1901 0 : ierr = pio_inq_varid(File, mdim_var_name, vardesc)
1902 0 : call cam_pio_handle_error(ierr, 'Error writing values for nonexistent mdimnames variable in write_hist_coord_vars')
1903 : ! Write out the values for mdim names
1904 0 : ierr = pio_put_var(File, vardesc, mdimnames)
1905 0 : call cam_pio_handle_error(ierr, 'Error writing values for mdimnames variable in write_hist_coord_vars')
1906 0 : deallocate(mdimnames)
1907 : end if
1908 :
1909 : ! Back to I/O or die trying
1910 9216 : call pio_seterrorhandling(File, PIO_INTERNAL_ERROR)
1911 :
1912 9216 : end subroutine write_hist_coord_vars
1913 :
1914 14733312 : subroutine lookup_hist_coord_indices(mdimnames, mdimindicies)
1915 : ! Dummy arguments
1916 : character(len=*), intent(in) :: mdimnames(:)
1917 : integer, intent(out) :: mdimindicies(:)
1918 :
1919 : ! Local variables
1920 : integer :: i, j
1921 : integer :: cnt
1922 : character(len=120) :: errormsg
1923 : character(len=16) :: name
1924 :
1925 :
1926 14733312 : cnt = size(mdimnames)
1927 23620608 : mdimindicies = -1
1928 :
1929 :
1930 23620608 : do j=1,cnt
1931 8887296 : name = mdimnames(j)
1932 59169792 : do i = 1, registeredmdims
1933 44436480 : if(name .eq. hist_coords(i)%name) then
1934 8887296 : mdimindicies(j)=i
1935 : end if
1936 : end do
1937 : end do
1938 23620608 : do j = 1, cnt
1939 23620608 : if(mdimindicies(j) < 0) then
1940 0 : do i = 1, registeredmdims
1941 0 : print *,__FILE__,__LINE__,i,hist_coords(i)%name
1942 : end do
1943 0 : write(errormsg,*) 'Name ',mdimnames(j),' is not a registered history coordinate'
1944 0 : call endrun(errormsg)
1945 : end if
1946 : end do
1947 :
1948 14733312 : end subroutine lookup_hist_coord_indices
1949 :
1950 : ! Find the vertical dimension (if present) in dimnames and return its size
1951 : ! (which is the number of levels). Return -1 if not found
1952 : ! If dimnames is not present, search all of the registered history coords
1953 14733312 : integer function hist_coord_find_levels(dimnames) result(levels)
1954 : ! Dummy argument
1955 : character(len=*), optional, intent(in) :: dimnames(:)
1956 :
1957 : ! Local variables
1958 : integer i, index, dimcnt
1959 :
1960 14733312 : levels = -1 ! Error return value
1961 :
1962 14733312 : if (present(dimnames)) then
1963 14733312 : dimcnt = size(dimnames)
1964 : else
1965 0 : dimcnt = registeredmdims
1966 : end if
1967 :
1968 14733312 : do i = 1, dimcnt
1969 8887296 : if (present(dimnames)) then
1970 8887296 : index = get_hist_coord_index(trim(dimnames(i)))
1971 8887296 : if (i < 0) then
1972 0 : call endrun('hist_coord_find_levels: '//trim(dimnames(i))//' is not a registred history coordinate')
1973 : end if
1974 : else
1975 : index = i ! Just cycle through all the registered mdims
1976 : end if
1977 :
1978 14733312 : if (hist_coords(index)%vertical_coord) then
1979 8887296 : levels = hist_coords(index)%dimsize
1980 8887296 : exit
1981 : end if
1982 : end do
1983 :
1984 14733312 : end function hist_coord_find_levels
1985 :
1986 : !#######################################################################
1987 :
1988 0 : character(len=max_hcoordname_len) function hist_dimension_name(size)
1989 : ! Given a specific size value, return the first registered dimension name which matches the size, if it exists
1990 : ! Otherwise the name returned is blank
1991 :
1992 : integer, intent(in) :: size
1993 :
1994 : integer :: i
1995 :
1996 0 : hist_dimension_name = ''
1997 :
1998 0 : do i=1, registeredmdims
1999 0 : if(size == hist_coords(i)%dimsize) then
2000 0 : hist_dimension_name = hist_coords(i)%name
2001 0 : exit
2002 : end if
2003 : end do
2004 :
2005 0 : end function hist_dimension_name
2006 :
2007 : !#######################################################################
2008 :
2009 0 : subroutine hist_dimension_values_r8(name, rvalues, istart, istop, found)
2010 : ! Given the name of a dimension, return its (real) values in <rvalues>
2011 : ! If <istart> and <istop> are present, they are the beginning and ending
2012 : ! indices of the dimension values to return in <rvalues>. By default,
2013 : ! the entire array is copied.
2014 : ! If <found> is passed, return .true. if <name> is a defined dimension
2015 : ! with real values.
2016 :
2017 : ! Dummy arguments
2018 : character(len=*), intent(in) :: name
2019 : real(r8), intent(out) :: rvalues(:)
2020 : integer, optional, intent(in) :: istart
2021 : integer, optional, intent(in) :: istop
2022 : logical, optional, intent(out) :: found
2023 : ! Local variables
2024 : integer :: indx, jndx, rndx
2025 : integer :: ibeg
2026 : integer :: iend
2027 : logical :: dim_ok
2028 : real(r8), parameter :: unset_r8 = huge(1.0_r8)
2029 : character(len=*), parameter :: subname = ': hist_dimension_values_r8'
2030 :
2031 0 : dim_ok = .false.
2032 0 : rvalues(:) = unset_r8
2033 :
2034 0 : do indx = 1, registeredmdims
2035 0 : if(trim(name) == trim(hist_coords(indx)%name)) then
2036 0 : dim_ok = associated(hist_coords(indx)%real_values)
2037 0 : if (dim_ok) then
2038 0 : if (present(istart)) then
2039 0 : ibeg = istart
2040 0 : if (ibeg < LBOUND(hist_coords(indx)%real_values, 1)) then
2041 0 : call endrun(subname//": istart is outside the bounds")
2042 : end if
2043 : else
2044 0 : ibeg = LBOUND(hist_coords(indx)%real_values, 1)
2045 : end if
2046 0 : if (present(istop)) then
2047 0 : iend = istop
2048 0 : if (iend > UBOUND(hist_coords(indx)%real_values, 1)) then
2049 0 : call endrun(subname//": istop is outside the bounds")
2050 : end if
2051 : else
2052 0 : iend = UBOUND(hist_coords(indx)%real_values, 1)
2053 : end if
2054 0 : if (SIZE(rvalues) < (iend - ibeg + 1)) then
2055 0 : call endrun(subname//": rvalues too small")
2056 : end if
2057 0 : rndx = 1
2058 0 : do jndx = ibeg, iend
2059 0 : rvalues(rndx) = hist_coords(indx)%real_values(jndx)
2060 0 : rndx = rndx + 1
2061 : end do
2062 : end if
2063 0 : exit
2064 : end if
2065 : end do
2066 0 : if (present(found)) then
2067 0 : found = dim_ok
2068 : end if
2069 :
2070 0 : end subroutine hist_dimension_values_r8
2071 :
2072 : !#######################################################################
2073 :
2074 0 : subroutine hist_dimension_values_int(name, ivalues, istart, istop, found)
2075 : ! Given the name of a dimension, return its (integer) values in <ivalues>
2076 : ! If <istart> and <istop> are present, they are the beginning and ending
2077 : ! indices of the dimension values to return in <ivalues>. By default,
2078 : ! the entire array is copied.
2079 : ! If <found> is passed, return .true. if <name> is a defined dimension
2080 : ! with integer values.
2081 :
2082 : ! Dummy arguments
2083 : character(len=*), intent(in) :: name
2084 : integer, intent(out) :: ivalues(:)
2085 : integer, optional, intent(in) :: istart
2086 : integer, optional, intent(in) :: istop
2087 : logical, optional, intent(out) :: found
2088 : ! Local variables
2089 : integer :: indx, jndx, rndx
2090 : integer :: ibeg
2091 : integer :: iend
2092 : logical :: dim_ok
2093 : integer, parameter :: unset_i = huge(1)
2094 : character(len=*), parameter :: subname = 'hist_dimension_values_int'
2095 :
2096 0 : dim_ok = .false.
2097 0 : ivalues(:) = unset_i
2098 :
2099 0 : do indx = 1, registeredmdims
2100 0 : if(trim(name) == trim(hist_coords(indx)%name)) then
2101 0 : dim_ok = associated(hist_coords(indx)%integer_values)
2102 0 : if (dim_ok) then
2103 0 : if (present(istart)) then
2104 0 : ibeg = istart
2105 0 : if (ibeg < LBOUND(hist_coords(indx)%integer_values, 1)) then
2106 0 : call endrun(subname//": istart is outside the bounds")
2107 : end if
2108 : else
2109 0 : ibeg = LBOUND(hist_coords(indx)%integer_values, 1)
2110 : end if
2111 0 : if (present(istop)) then
2112 0 : iend = istop
2113 0 : if (iend > UBOUND(hist_coords(indx)%integer_values, 1)) then
2114 0 : call endrun(subname//": istop is outside the bounds")
2115 : end if
2116 : else
2117 0 : iend = UBOUND(hist_coords(indx)%integer_values, 1)
2118 : end if
2119 0 : if (SIZE(ivalues) < (iend - ibeg + 1)) then
2120 0 : call endrun(subname//": ivalues too small")
2121 : end if
2122 0 : rndx = 1
2123 0 : do jndx = ibeg, iend
2124 0 : ivalues(rndx) = hist_coords(indx)%integer_values(jndx)
2125 0 : rndx = rndx + 1
2126 : end do
2127 : end if
2128 0 : exit
2129 : end if
2130 : end do
2131 0 : if (present(found)) then
2132 0 : found = dim_ok
2133 : end if
2134 :
2135 0 : end subroutine hist_dimension_values_int
2136 :
2137 : !#######################################################################
2138 :
2139 0 : end module cam_history_support
|