14 use netcdf,
only : nf90_noerr, nf90_global, nf90_nowrite, nf90_inquire_attribute, nf90_open, nf90_strerror, &
15 nf90_inq_dimid, nf90_inquire_dimension, nf90_inq_varid, nf90_get_var, nf90_inquire, nf90_close, nf90_get_att
23 character(len=*),
parameter ::
time_key =
"time",& !< Corresponding NetCDF data time key
55 type(model_state_type),
target,
intent(inout) :: current_state
57 integer :: ncid, time_dim, z_dim, z_half_dim, x_dim, x_half_dim
58 real,
dimension(:),
allocatable :: time, x, z, x_half, z_half
59 real,
dimension(:,:,:),
allocatable :: u, w, v
61 if (options_has_key(current_state%options_database,
"restart"))
then
62 call log_master_log(log_debug,
"Ignoring KiD reader as restart from checkpoint file selected")
66 configuration_file=options_get_string(current_state%options_database,
"kid_configuration_file")
67 current_state%rhobous=options_get_real(current_state%options_database,
"rhobous")
68 current_state%thref0=options_get_real(current_state%options_database,
"thref0")
69 current_state%surface_pressure=options_get_real(current_state%options_database,
"surface_pressure")
70 flood_q=options_get_logical(current_state%options_database,
"flood_q")
71 float_q=options_get_logical(current_state%options_database,
"float_q")
72 clone_to_3d=options_get_logical(current_state%options_database,
"clone_to_3d")
73 rotate_xy=options_get_logical(current_state%options_database,
"rotate_xy")
76 call options_get_integer_array(current_state%options_database,
"q_coordinates_x",
q_coordinates_x)
77 call options_get_integer_array(current_state%options_database,
"q_coordinates_y",
q_coordinates_y)
78 call options_get_integer_array(current_state%options_database,
"q_coordinates_z",
q_coordinates_z)
79 call options_get_integer_array(current_state%options_database,
"q_coordinates_value",
q_coordinates_value)
83 if (log_get_logging_level() .ge. log_debug)
call read_global_attributes(ncid, current_state%parallel%my_rank)
84 call read_dimensions(ncid, time_dim, z_dim, z_half_dim, x_dim, x_half_dim)
85 call read_variables(ncid, time_dim, z_dim, z_half_dim, x_dim, x_half_dim, time, x, z, x_half, z_half, u, w, v)
88 call create_grid(current_state%global_grid, z_half, x_half, z_half_dim, x_half_dim, current_state%parallel%my_rank)
102 current_state%initialised=.true.
103 call log_master_log(log_info,
"Initialised configuration from KiD model file `"//trim(
configuration_file)//
"`")
107 type(model_state_type),
intent(inout) :: current_state
109 integer :: x_size, y_size, z_size, i
111 z_size = current_state%local_grid%size(z_index) + current_state%local_grid%halo_size(z_index) * 2
112 y_size = current_state%local_grid%size(y_index) + current_state%local_grid%halo_size(y_index) * 2
113 x_size = current_state%local_grid%size(x_index) + current_state%local_grid%halo_size(x_index) * 2
115 if (
flood_q) current_state%number_q_fields=current_state%number_q_fields+1
116 if (
float_q) current_state%number_q_fields=current_state%number_q_fields+1
117 allocate(current_state%q(current_state%number_q_fields), current_state%zq(current_state%number_q_fields), &
118 current_state%sq(current_state%number_q_fields))
119 do i=1,current_state%number_q_fields
121 if (
flood_q .and. i == 1) current_state%q(i)%data(:,:,:) = 1.
132 type(model_state_type),
intent(inout) :: current_state
133 type(prognostic_field_type),
intent(inout) :: q_field
135 integer :: i, x_local, y_local
143 x_local=(
q_coordinates_x(i) - (current_state%local_grid%start(x_index)-1)) + current_state%local_grid%halo_size(x_index)
144 y_local=(
q_coordinates_y(i) - (current_state%local_grid%start(y_index)-1)) + current_state%local_grid%halo_size(y_index)
146 q_field%data(:,current_state%local_grid%local_domain_start_index(y_index):&
147 current_state%local_grid%local_domain_end_index(y_index),current_state%local_grid%local_domain_start_index(x_index):&
150 q_field%data(:,current_state%local_grid%local_domain_start_index(y_index):&
153 q_field%data(:,y_local,current_state%local_grid%local_domain_start_index(x_index):&
158 q_field%data(
q_coordinates_z(i),current_state%local_grid%local_domain_start_index(y_index):&
159 current_state%local_grid%local_domain_end_index(y_index),current_state%local_grid%local_domain_start_index(x_index):&
162 q_field%data(
q_coordinates_z(i),current_state%local_grid%local_domain_start_index(y_index):&
165 q_field%data(
q_coordinates_z(i),y_local,current_state%local_grid%local_domain_start_index(x_index):&
175 type(model_state_type),
intent(inout) :: current_state
176 integer,
intent(in) :: q_id, x_size, y_size, z_size
178 allocate(current_state%q(q_id)%data(z_size, y_size, x_size), current_state%zq(q_id)%data(z_size, y_size, x_size), &
179 current_state%sq(q_id)%data(z_size, y_size, x_size))
180 current_state%q(q_id)%data(:,:,:) = 0.
181 current_state%zq(q_id)%data(:,:,:) = 0.
182 current_state%sq(q_id)%data(:,:,:) = 0.
183 current_state%q(q_id)%active=.true.
184 current_state%zq(q_id)%active=.true.
185 current_state%sq(q_id)%active=.true.
192 type(model_state_type),
intent(inout) :: current_state
194 if (
associated(current_state%parallel%decomposition_procedure))
then
195 call current_state%parallel%decomposition_procedure(current_state)
197 call log_master_log(log_error,
"No decomposition specified")
204 type(model_state_type),
intent(inout) :: current_state
206 integer :: x_size, y_size, z_size
208 z_size = current_state%local_grid%size(z_index) + current_state%local_grid%halo_size(z_index) * 2
209 y_size = current_state%local_grid%size(y_index) + current_state%local_grid%halo_size(y_index) * 2
210 x_size = current_state%local_grid%size(x_index) + current_state%local_grid%halo_size(x_index) * 2
213 allocate(current_state%zu%data(z_size, y_size, x_size))
214 allocate(current_state%su%data(z_size, y_size, x_size))
215 allocate(current_state%savu%data(z_size, y_size, x_size))
216 current_state%zu%data(:,:,:)= 0.0
219 allocate(current_state%zv%data(z_size, y_size, x_size))
220 allocate(current_state%sv%data(z_size, y_size, x_size))
221 allocate(current_state%savv%data(z_size, y_size, x_size))
222 current_state%zv%data(:,:,:)= 0.0
225 allocate(current_state%zw%data(z_size, y_size, x_size))
226 allocate(current_state%sw%data(z_size, y_size, x_size))
227 allocate(current_state%savw%data(z_size, y_size, x_size))
228 current_state%zw%data(:,:,:)= 0.0
230 if (current_state%th%active)
then
231 allocate(current_state%zth%data(z_size, y_size, x_size))
232 allocate(current_state%sth%data(z_size, y_size, x_size))
233 current_state%zth%data(:,:,:)= 0.0
245 type(local_grid_type),
intent(inout) :: local_grid
246 type(prognostic_field_type),
intent(inout) :: field
247 integer,
intent(in) :: z_grid, y_grid, x_grid
248 real,
dimension(:,:,:),
allocatable,
intent(in) :: data
250 integer :: i, j, k, preMulSizeY, preMulSizeX
252 field%grid(z_index) = z_grid
253 field%grid(y_index) = y_grid
254 field%grid(x_index) = x_grid
255 field%active = .true.
257 allocate(field%data(local_grid%size(z_index) + local_grid%halo_size(z_index) * 2, local_grid%size(y_index) + &
258 local_grid%halo_size(y_index) * 2, local_grid%size(x_index) + local_grid%halo_size(x_index) * 2))
266 do k=local_grid%start(z_index),local_grid%end(z_index)
267 field%data(local_grid%halo_size(z_index)+(k-((local_grid%start(z_index)-1))), local_grid%halo_size(y_index)+&
268 (j-((ceiling(local_grid%start(y_index)/ real(
domain_multiplication))-1))), local_grid%halo_size(x_index)+&
270 real(data(1, k, merge(j, i, rotate_xy)), kind=default_precision)
279 if (local_grid%active(y_index))
then
281 field%data(:, local_grid%halo_size(y_index) + i*premulsizey +1 : local_grid%halo_size(y_index) + (i+1)*premulsizey, &
282 local_grid%halo_size(x_index)+1: local_grid%halo_size(x_index)+premulsizex) = &
283 field%data(:, local_grid%halo_size(y_index)+1 : local_grid%halo_size(y_index) + premulsizey, &
284 local_grid%halo_size(x_index)+1 : local_grid%halo_size(x_index)+premulsizex)
287 if (local_grid%active(x_index))
then
289 field%data(:, local_grid%local_domain_start_index(y_index) : local_grid%local_domain_end_index(y_index), &
290 local_grid%halo_size(x_index) + i*premulsizex + 1 : local_grid%halo_size(x_index) + (i+1)*premulsizex) = &
291 field%data(:, local_grid%local_domain_start_index(y_index) : local_grid%local_domain_end_index(y_index), &
292 local_grid%halo_size(x_index)+1 : local_grid%halo_size(x_index)+premulsizex)
304 subroutine create_grid(specific_grid, z, x, z_dim, x_dim, my_rank)
305 type(global_grid_type),
intent(inout) :: specific_grid
306 integer,
intent(in) :: z_dim, x_dim, my_rank
307 real,
dimension(:),
intent(in) :: z, x
309 specific_grid%bottom(z_index) = int(z(1))
310 specific_grid%bottom(merge(y_index, x_index,
rotate_xy)) = int(x(1))
311 if (
clone_to_3d) specific_grid%bottom(y_index) = int(x(1))
313 specific_grid%top(z_index) = int(z(z_dim))
317 specific_grid%resolution(z_index) = int(z(2) - z(1))
318 specific_grid%resolution(merge(y_index, x_index,
rotate_xy)) = int(x(2) - x(1))
319 if (
clone_to_3d) specific_grid%resolution(y_index) = specific_grid%resolution(x_index)
321 specific_grid%size(z_index) = z_dim
325 specific_grid%active(z_index) = .true.
326 specific_grid%active(merge(y_index, x_index,
rotate_xy)) = .true.
327 if (
clone_to_3d) specific_grid%active(y_index) = .true.
330 if (.not. specific_grid%active(x_index))
call log_master_log(log_error, &
331 "Model compiled with X active but inactive in configuration")
333 if (specific_grid%active(x_index))
call log_master_log(log_error, &
334 "Model compiled with X inactive but active in configuration")
337 if (.not. specific_grid%active(y_index))
call log_master_log(log_error, &
338 "Model compiled with Y active but inactive in configuration")
340 if (specific_grid%active(y_index))
call log_master_log(log_error, &
341 "Model compiled with Y inactive but active in configuration")
344 if (.not. specific_grid%active(z_index))
call log_master_log(log_error, &
345 "Model compiled with Z active but inactive in configuration")
347 if (specific_grid%active(z_index))
call log_master_log(log_error, &
348 "Model compiled with Z inactive but active in configuration")
350 specific_grid%dimensions = merge(3, 2,
clone_to_3d)
359 type(model_state_type),
intent(inout) :: current_state
360 integer,
intent(in) :: z_size
361 real,
dimension(:),
intent(in) :: z
365 allocate(current_state%global_grid%configuration%vertical%kgd(z_size), &
366 current_state%global_grid%configuration%vertical%hgd(z_size))
369 current_state%global_grid%configuration%vertical%kgd(i) = i
370 current_state%global_grid%configuration%vertical%hgd(i) = real(z(i))
379 integer,
intent(in) :: ncid
381 integer :: ndims_in, nvars_in, ngatts_in, unlimdimid_in
383 call check_status(nf90_inquire(ncid, ndims_in, nvars_in, ngatts_in, unlimdimid_in))
384 if (ndims_in /= 5)
call log_log(log_error,
"NetCDF KiD number of model dimensions must equal 5")
385 if (nvars_in /= 7)
call log_log(log_error,
"NetCDF KiD number of model variables must equal 5")
386 if (ngatts_in .le. 0)
call log_log(log_error,
"NetCDF KiD global attributes must be specified")
387 if (unlimdimid_in .gt. 0)
call log_log(log_error,
"NetCDF KiD model number of unlimited dimensions must be 0")
404 subroutine read_variables(ncid, time_dim, z_dim, z_half_dim, x_dim, x_half_dim, time, x, z, x_half, z_half, u, w, v)
405 integer,
intent(in) :: ncid, time_dim, z_dim, z_half_dim, x_dim, x_half_dim
406 real,
dimension(:),
allocatable,
intent(inout) :: time, x, z, x_half, z_half
407 real,
dimension(:,:,:),
allocatable,
intent(inout) :: u, w, v
409 allocate(time(time_dim))
412 allocate(z_half(z_half_dim))
413 allocate(x_half(x_half_dim))
416 allocate(v(time_dim, z_dim, x_half_dim))
418 allocate(u(time_dim, z_dim, x_half_dim))
420 allocate(w(time_dim, z_half_dim, x_half_dim))
435 if (.not.
allocated(v))
then
436 allocate(v(time_dim, z_dim, x_half_dim))
439 if (.not.
allocated(u))
then
440 allocate(u(time_dim, z_dim, x_half_dim))
452 integer,
intent(in) :: ncid
453 character(len=*),
intent(in) :: key
454 real,
dimension(:),
intent(inout),
optional :: data1d
455 real,
dimension(:,:,:),
intent(inout),
optional :: data3d
457 integer :: variable_id
458 real,
dimension(:,:,:),
allocatable :: sdata
460 call check_status(nf90_inq_varid(ncid, key, variable_id))
462 if (.not.
present(data1d) .and. .not.
present(data3d))
return
464 if (
present(data1d))
then
465 call check_status(nf90_get_var(ncid, variable_id, data1d))
468 allocate(sdata(
size(data3d,1),
size(data3d,3),
size(data3d,2)))
469 call check_status(nf90_get_var(ncid, variable_id, sdata))
470 data3d(:,:,:)=reshape(sdata(:,:,:),(/
size(data3d,1),
size(data3d,2),
size(data3d,3)/))
483 integer,
intent(in) :: ncid
484 integer,
intent(out) :: time_dim, z_dim, z_half_dim, x_dim, x_half_dim
486 integer :: time_dimid, z_dimid, z_half_dimid, x_dimid, x_half_dimid
494 call check_status(nf90_inquire_dimension(ncid, time_dimid, len=time_dim))
495 call check_status(nf90_inquire_dimension(ncid, z_dimid, len=z_dim))
496 call check_status(nf90_inquire_dimension(ncid, z_half_dimid, len=z_half_dim))
497 call check_status(nf90_inquire_dimension(ncid, x_dimid, len=x_dim))
498 call check_status(nf90_inquire_dimension(ncid, x_half_dimid, len=x_half_dim))
504 integer,
intent(in) :: ncid, pid
506 character(len=:),
allocatable :: attributeValue
509 call log_master_log(log_debug,
"KiD file title: "//attributevalue)
510 deallocate(attributevalue)
513 call log_master_log(log_debug,
"KiD file created: "//attributevalue)
514 deallocate(attributevalue)
522 integer,
intent(in) :: ncid
523 character(len=*),
intent(in) :: key
528 call check_status(nf90_inquire_attribute(ncid, nf90_global, key, len = length))
536 integer,
intent(in) :: status
538 if (status /= nf90_noerr)
then
539 call log_log(log_error,
"NetCDF returned error code of "//trim(nf90_strerror(status)))