MONC
Functions/Subroutines | Variables
kidreader_mod Module Reference

Component to set up the model based upon a KiD model configuration. More...

Functions/Subroutines

type(component_descriptor_type) function, public kidreader_get_descriptor ()
 Provides the descriptor back to the caller and is used in component registration. More...
 
subroutine initialise_callback (current_state)
 Initialisation hook which will parse the configuration NetCDF file and set up the model based upon this. More...
 
subroutine set_up_q_fields (current_state)
 
subroutine populate_q_tracer (current_state, q_field)
 Populates the Q tracer field based upon the configuration that has been read in from the simulation file. More...
 
subroutine initialise_single_q_field (current_state, q_id, z_size, y_size, x_size)
 
subroutine decompose_grid (current_state)
 Calls the decomposition procedure to decompose the grid and determine neighbouring processes If no decomposition procedure is specified then this results in an error. More...
 
subroutine initalise_source_and_z_fields (current_state)
 Based upon the local grid this will initialise the Source, Z and SAV fields for each prognostic. More...
 
subroutine initialise_velocity_field (local_grid, field, z_grid, y_grid, x_grid, data)
 Will initialise a velocity field with the loaded data. More...
 
subroutine create_grid (specific_grid, z, x, z_dim, x_dim, my_rank)
 Creates a specific grid based upon the data read from the KiD model NetCDF file. More...
 
subroutine define_vertical_levels (current_state, z, z_size)
 Defines the vertical levels of the grid. This is both the grid points and corresponding height for each point in metres. More...
 
subroutine check_kinematics_file (ncid)
 Checks that the kinematics file that has been loaded is consistent with what we expect. More...
 
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)
 Reads the variables from the NetCDF KiD model file. More...
 
subroutine read_single_variable (ncid, key, data1d, data3d)
 Reads a single variable out of a NetCDF file. More...
 
subroutine read_dimensions (ncid, time_dim, z_dim, z_half_dim, x_dim, x_half_dim)
 Reads the dimensions from the NetCDF file. More...
 
subroutine read_global_attributes (ncid, pid)
 Will read the global attributes of the NetCDF KiD model dump and log_log them to debug. More...
 
character(len=:) function, allocatable, target read_specific_global_attribute (ncid, key)
 Will read a specific global NetCDF attribute. More...
 
subroutine check_status (status)
 Will check a NetCDF status and write to log_log error any decoded statuses. More...
 

Variables

integer, parameter number_q_coords = 100
 Number of Q field value coords that can be specified. More...
 
character(len= *), parameter time_key = "time"
 Corresponding NetCDF data time key. More...
 
character(len= *), parameter z_key = "z"
 Corresponding NetCDF data z (primal grid) key. More...
 
character(len= *), parameter z_half_key = "z_half"
 Corresponding NetCDF data z half (dual grid) key. More...
 
character(len= *), parameter x_key = "x"
 Corresponding NetCDF data x (primal grid) key. More...
 
character(len= *), parameter x_half_key = "x_half"
 Corresponding NetCDF data x half (primal grid) key. More...
 
character(len= *), parameter u_key ="u"
 Corresponding NetCDF data u flow field key. More...
 
character(len= *), parameter w_key ="w"
 Corresponding NetCDF data w flow field key. More...
 
logical flood_q = .false.
 
logical float_q = .false.
 
logical clone_to_3d = .false.
 
logical rotate_xy =.false.
 
integer domain_multiplication =1
 The Q field coordinates configured by the user. More...
 
integer, dimension(number_q_coords), save q_coordinates_x
 
integer, dimension(number_q_coords), save q_coordinates_y
 
integer, dimension(number_q_coords), save q_coordinates_z
 
integer, dimension(number_q_coords), save q_coordinates_value
 
character(len=string_length) configuration_file
 NetCDF model file to load. More...
 

Detailed Description

Component to set up the model based upon a KiD model configuration.

These data files are in NetCDF format

Function/Subroutine Documentation

◆ check_kinematics_file()

subroutine kidreader_mod::check_kinematics_file ( integer, intent(in)  ncid)
private

Checks that the kinematics file that has been loaded is consistent with what we expect.

This checks the file against a number of limitations of the current KiD reading functionality to ensure that this component can parse it correctly

Definition at line 378 of file kidreader.F90.

379  integer, intent(in) :: ncid
380 
381  integer :: ndims_in, nvars_in, ngatts_in, unlimdimid_in
382 
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")
Here is the call graph for this function:
Here is the caller graph for this function:

◆ check_status()

subroutine kidreader_mod::check_status ( integer, intent(in)  status)
private

Will check a NetCDF status and write to log_log error any decoded statuses.

Parameters
statusThe NetCDF status flag

Definition at line 535 of file kidreader.F90.

536  integer, intent(in) :: status
537 
538  if (status /= nf90_noerr) then
539  call log_log(log_error, "NetCDF returned error code of "//trim(nf90_strerror(status)))
540  end if
Here is the caller graph for this function:

◆ create_grid()

subroutine kidreader_mod::create_grid ( type(global_grid_type), intent(inout)  specific_grid,
real, dimension(:), intent(in)  z,
real, dimension(:), intent(in)  x,
integer, intent(in)  z_dim,
integer, intent(in)  x_dim,
integer, intent(in)  my_rank 
)
private

Creates a specific grid based upon the data read from the KiD model NetCDF file.

Parameters
specific_gridThe grid to create
zArray of grid points in the z dimension
xArray of grid points in the x dimension
z_dimThe number of points in the z dimension
x_dimThe number of points in the x dimension

Definition at line 304 of file kidreader.F90.

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
308 
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))
312 
313  specific_grid%top(z_index) = int(z(z_dim))
314  specific_grid%top(merge(y_index, x_index, rotate_xy)) = int(x(x_dim)) * domain_multiplication
315  if (clone_to_3d) specific_grid%top(y_index) = int(x(x_dim)) * domain_multiplication
316 
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)
320 
321  specific_grid%size(z_index) = z_dim
322  specific_grid%size(merge(y_index, x_index, rotate_xy)) = x_dim * domain_multiplication
323  if (clone_to_3d) specific_grid%size(y_index) = x_dim * domain_multiplication
324 
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.
328 
329 #ifdef U_ACTIVE
330  if (.not. specific_grid%active(x_index)) call log_master_log(log_error, &
331  "Model compiled with X active but inactive in configuration")
332 #else
333  if (specific_grid%active(x_index)) call log_master_log(log_error, &
334  "Model compiled with X inactive but active in configuration")
335 #endif
336 #ifdef V_ACTIVE
337  if (.not. specific_grid%active(y_index)) call log_master_log(log_error, &
338  "Model compiled with Y active but inactive in configuration")
339 #else
340  if (specific_grid%active(y_index)) call log_master_log(log_error, &
341  "Model compiled with Y inactive but active in configuration")
342 #endif
343 #ifdef W_ACTIVE
344  if (.not. specific_grid%active(z_index)) call log_master_log(log_error, &
345  "Model compiled with Z active but inactive in configuration")
346 #else
347  if (specific_grid%active(z_index)) call log_master_log(log_error, &
348  "Model compiled with Z inactive but active in configuration")
349 #endif
350  specific_grid%dimensions = merge(3, 2, clone_to_3d)
Here is the caller graph for this function:

◆ decompose_grid()

subroutine kidreader_mod::decompose_grid ( type(model_state_type), intent(inout)  current_state)
private

Calls the decomposition procedure to decompose the grid and determine neighbouring processes If no decomposition procedure is specified then this results in an error.

Parameters
current_stateThe current model state_mod

Definition at line 191 of file kidreader.F90.

192  type(model_state_type), intent(inout) :: current_state
193 
194  if (associated(current_state%parallel%decomposition_procedure)) then
195  call current_state%parallel%decomposition_procedure(current_state)
196  else
197  call log_master_log(log_error, "No decomposition specified")
198  end if
Here is the caller graph for this function:

◆ define_vertical_levels()

subroutine kidreader_mod::define_vertical_levels ( type(model_state_type), intent(inout)  current_state,
real, dimension(:), intent(in)  z,
integer, intent(in)  z_size 
)
private

Defines the vertical levels of the grid. This is both the grid points and corresponding height for each point in metres.

Parameters
current_stateThe current model state_mod
zArray of grid heights in metres
z_sizeNumber of grid points

Definition at line 358 of file kidreader.F90.

359  type(model_state_type), intent(inout) :: current_state
360  integer, intent(in) :: z_size
361  real, dimension(:), intent(in) :: z
362 
363  integer :: i
364 
365  allocate(current_state%global_grid%configuration%vertical%kgd(z_size), &
366  current_state%global_grid%configuration%vertical%hgd(z_size))
367 
368  do i=1,z_size
369  current_state%global_grid%configuration%vertical%kgd(i) = i
370  current_state%global_grid%configuration%vertical%hgd(i) = real(z(i))
371  end do
Here is the caller graph for this function:

◆ initalise_source_and_z_fields()

subroutine kidreader_mod::initalise_source_and_z_fields ( type(model_state_type), intent(inout)  current_state)
private

Based upon the local grid this will initialise the Source, Z and SAV fields for each prognostic.

Parameters
current_state

Definition at line 203 of file kidreader.F90.

204  type(model_state_type), intent(inout) :: current_state
205 
206  integer :: x_size, y_size, z_size
207 
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
211 
212 #ifdef U_ACTIVE
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
217 #endif
218 #ifdef V_ACTIVE
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
223 #endif
224 #ifdef W_ACTIVE
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
229 #endif
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
234  end if
Here is the caller graph for this function:

◆ initialise_callback()

subroutine kidreader_mod::initialise_callback ( type(model_state_type), intent(inout), target  current_state)
private

Initialisation hook which will parse the configuration NetCDF file and set up the model based upon this.

Parameters
current_stateThe current model state_mod

Definition at line 54 of file kidreader.F90.

55  type(model_state_type), target, intent(inout) :: current_state
56 
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
60 
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")
63  return
64  end if
65 
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")
74  domain_multiplication=options_get_integer(current_state%options_database, "domain_multiplication")
75 
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)
80 
81  call check_status(nf90_open(path = trim(configuration_file), mode = nf90_nowrite, ncid = ncid))
82  call check_kinematics_file(ncid)
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)
86  call check_status(nf90_close(ncid))
87 
88  call create_grid(current_state%global_grid, z_half, x_half, z_half_dim, x_half_dim, current_state%parallel%my_rank)
89  call define_vertical_levels(current_state, z_half, z_half_dim)
90 
91  call decompose_grid(current_state)
92 
93  if (rotate_xy) then
94  call initialise_velocity_field(current_state%local_grid, current_state%v, dual_grid, primal_grid, dual_grid, v)
95  else
96  call initialise_velocity_field(current_state%local_grid, current_state%u, dual_grid, dual_grid, primal_grid, u)
97  end if
98  call initialise_velocity_field(current_state%local_grid, current_state%w, primal_grid, dual_grid, dual_grid, w)
99  if (clone_to_3d) call initialise_velocity_field(current_state%local_grid, current_state%v, dual_grid, primal_grid, dual_grid, v)
100  call initalise_source_and_z_fields(current_state)
101  call set_up_q_fields(current_state)
102  current_state%initialised=.true.
103  call log_master_log(log_info, "Initialised configuration from KiD model file `"//trim(configuration_file)//"`")
Here is the call graph for this function:
Here is the caller graph for this function:

◆ initialise_single_q_field()

subroutine kidreader_mod::initialise_single_q_field ( type(model_state_type), intent(inout)  current_state,
integer, intent(in)  q_id,
integer, intent(in)  z_size,
integer, intent(in)  y_size,
integer, intent(in)  x_size 
)
private

Definition at line 174 of file kidreader.F90.

175  type(model_state_type), intent(inout) :: current_state
176  integer, intent(in) :: q_id, x_size, y_size, z_size
177 
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.
Here is the caller graph for this function:

◆ initialise_velocity_field()

subroutine kidreader_mod::initialise_velocity_field ( type(local_grid_type), intent(inout)  local_grid,
type(prognostic_field_type), intent(inout)  field,
integer, intent(in)  z_grid,
integer, intent(in)  y_grid,
integer, intent(in)  x_grid,
real, dimension(:,:,:), intent(in), allocatable  data 
)
private

Will initialise a velocity field with the loaded data.

Parameters
fieldThe velocity field to initialise
z_dimThe size of the z dimension
y_dimThe size of the y dimension
x_dimThe size of the x dimension
gridThe applicable grid to place this field upon
dataThe NetCDF KiD model data to load into the velocity field

Definition at line 244 of file kidreader.F90.

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
249 
250  integer :: i, j, k, preMulSizeY, preMulSizeX
251 
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.
256 
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))
259  field%data=0.0
260 
261  ! Divisions here are for the multiplication - we just fill up the original size and then duplicate this across dimensions
262  do i=ceiling(local_grid%start(x_index)/real(domain_multiplication)),&
263  ceiling(local_grid%end(x_index)/real(domain_multiplication))
264  do j=ceiling(local_grid%start(y_index)/real(domain_multiplication)),&
265  ceiling(local_grid%end(y_index)/real(domain_multiplication))
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)+&
269  (i- ((ceiling(local_grid%start(x_index)/real(domain_multiplication))-1)))) = &
270  real(data(1, k, merge(j, i, rotate_xy)), kind=default_precision)
271  end do
272  end do
273  end do
274 
275  if (domain_multiplication .ge. 2) then
276  ! If there is a domain multiplication then duplicate data across dimensions
277  premulsizey = ceiling(local_grid%size(y_index) / real(domain_multiplication))
278  premulsizex = ceiling(local_grid%size(x_index) / real(domain_multiplication))
279  if (local_grid%active(y_index)) then
280  do i=1,domain_multiplication-1
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)
285  end do
286  end if
287  if (local_grid%active(x_index)) then
288  do i=1,domain_multiplication-1
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)
293  end do
294  end if
295  end if
Here is the caller graph for this function:

◆ kidreader_get_descriptor()

type(component_descriptor_type) function, public kidreader_mod::kidreader_get_descriptor

Provides the descriptor back to the caller and is used in component registration.

Returns
The KidReader component descriptor

Definition at line 45 of file kidreader.F90.

46  kidreader_get_descriptor%name="kidreader"
47  kidreader_get_descriptor%version=0.1
48  kidreader_get_descriptor%initialisation=>initialise_callback
Here is the call graph for this function:

◆ populate_q_tracer()

subroutine kidreader_mod::populate_q_tracer ( type(model_state_type), intent(inout)  current_state,
type(prognostic_field_type), intent(inout)  q_field 
)
private

Populates the Q tracer field based upon the configuration that has been read in from the simulation file.

Parameters
current_stateThe current model state_mod
q_fieldThe qfield that we are going to fill in

Definition at line 131 of file kidreader.F90.

132  type(model_state_type), intent(inout) :: current_state
133  type(prognostic_field_type), intent(inout) :: q_field
134 
135  integer :: i, x_local, y_local
136 
137  do i=1,number_q_coords
138  if (q_coordinates_x(i) .ne. -1 .and. q_coordinates_y(i) .ne. -1) then
139  if (q_coordinates_x(i) .ne. 0 .and. .not. (q_coordinates_x(i) .ge. current_state%local_grid%start(x_index) .and. &
140  q_coordinates_x(i) .le. current_state%local_grid%end(x_index))) cycle
141  if (q_coordinates_y(i) .ne. 0 .and. .not. (q_coordinates_y(i) .ge. current_state%local_grid%start(y_index) .and. &
142  q_coordinates_y(i) .le. current_state%local_grid%end(y_index))) cycle
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)
145  if (q_coordinates_z(i) == 0 .and. q_coordinates_y(i) == 0 .and. q_coordinates_x(i) == 0) then
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):&
148  current_state%local_grid%local_domain_end_index(x_index)) = q_coordinates_value(i)
149  else if (q_coordinates_z(i) == 0 .and. q_coordinates_y(i) == 0 .and. q_coordinates_x(i) .ne. 0) then
150  q_field%data(:,current_state%local_grid%local_domain_start_index(y_index):&
151  current_state%local_grid%local_domain_end_index(y_index),x_local) = q_coordinates_value(i)
152  else if (q_coordinates_z(i) == 0 .and. q_coordinates_y(i) .ne. 0 .and. q_coordinates_x(i) == 0) then
153  q_field%data(:,y_local,current_state%local_grid%local_domain_start_index(x_index):&
154  current_state%local_grid%local_domain_end_index(x_index)) = q_coordinates_value(i)
155  else if (q_coordinates_z(i) == 0 .and. q_coordinates_y(i) .ne. 0 .and. q_coordinates_x(i) .ne. 0) then
156  q_field%data(:,y_local,x_local) = q_coordinates_value(i)
157  else if (q_coordinates_z(i) .ne. 0 .and. q_coordinates_y(i) == 0 .and. q_coordinates_x(i) == 0) then
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):&
160  current_state%local_grid%local_domain_end_index(x_index)) = q_coordinates_value(i)
161  else if (q_coordinates_z(i) .ne. 0 .and. q_coordinates_y(i) == 0 .and. q_coordinates_x(i) .ne. 0) then
162  q_field%data(q_coordinates_z(i),current_state%local_grid%local_domain_start_index(y_index):&
163  current_state%local_grid%local_domain_end_index(y_index),x_local) = q_coordinates_value(i)
164  else if (q_coordinates_z(i) .ne. 0 .and. q_coordinates_y(i) .ne. 0 .and. q_coordinates_x(i) == 0) then
165  q_field%data(q_coordinates_z(i),y_local,current_state%local_grid%local_domain_start_index(x_index):&
166  current_state%local_grid%local_domain_end_index(x_index)) = q_coordinates_value(i)
167  else if (q_coordinates_z(i) .ne. 0 .and. q_coordinates_y(i) .ne. 0 .and. q_coordinates_x(i) .ne. 0) then
168  q_field%data(q_coordinates_z(i),y_local,x_local) = q_coordinates_value(i)
169  end if
170  end if
171  end do
Here is the caller graph for this function:

◆ read_dimensions()

subroutine kidreader_mod::read_dimensions ( integer, intent(in)  ncid,
integer, intent(out)  time_dim,
integer, intent(out)  z_dim,
integer, intent(out)  z_half_dim,
integer, intent(out)  x_dim,
integer, intent(out)  x_half_dim 
)
private

Reads the dimensions from the NetCDF file.

Parameters
ncidThe NetCDF file id
time_dimNumber of elements in the time dimension
z_dimNumber of elements in the z dimension of the primal grid
z_half_dimNumber of elements in the z dimension of the dual grid
x_dimNumber of elements in the x dimension of the primal grid
x_half_dimNumber of elements in the x dimension of the dual grid

Definition at line 482 of file kidreader.F90.

483  integer, intent(in) :: ncid
484  integer, intent(out) :: time_dim, z_dim, z_half_dim, x_dim, x_half_dim
485 
486  integer :: time_dimid, z_dimid, z_half_dimid, x_dimid, x_half_dimid
487 
488  call check_status(nf90_inq_dimid(ncid, time_key, time_dimid))
489  call check_status(nf90_inq_dimid(ncid, z_key, z_dimid))
490  call check_status(nf90_inq_dimid(ncid, z_half_key, z_half_dimid))
491  call check_status(nf90_inq_dimid(ncid, x_key, x_dimid))
492  call check_status(nf90_inq_dimid(ncid, x_half_key, x_half_dimid))
493 
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))
Here is the call graph for this function:
Here is the caller graph for this function:

◆ read_global_attributes()

subroutine kidreader_mod::read_global_attributes ( integer, intent(in)  ncid,
integer, intent(in)  pid 
)
private

Will read the global attributes of the NetCDF KiD model dump and log_log them to debug.

Parameters
ncidThe NetCDF file id

Definition at line 503 of file kidreader.F90.

504  integer, intent(in) :: ncid, pid
505 
506  character(len=:), allocatable :: attributeValue
507 
508  attributevalue=read_specific_global_attribute(ncid, "title")
509  call log_master_log(log_debug, "KiD file title: "//attributevalue)
510  deallocate(attributevalue)
511 
512  attributevalue=read_specific_global_attribute(ncid, "creation date")
513  call log_master_log(log_debug, "KiD file created: "//attributevalue)
514  deallocate(attributevalue)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ read_single_variable()

subroutine kidreader_mod::read_single_variable ( integer, intent(in)  ncid,
character(len=*), intent(in)  key,
real, dimension(:), intent(inout), optional  data1d,
real, dimension(:,:,:), intent(inout), optional  data3d 
)
private

Reads a single variable out of a NetCDF file.

Parameters
ncidThe NetCDF file id
keyThe variable key (name) to access
data1dOptional one dimensional data to read into
data3dOptional three dimensional data to read into

Definition at line 451 of file kidreader.F90.

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
456 
457  integer :: variable_id
458  real, dimension(:,:,:), allocatable :: sdata
459 
460  call check_status(nf90_inq_varid(ncid, key, variable_id))
461 
462  if (.not. present(data1d) .and. .not. present(data3d)) return
463 
464  if (present(data1d)) then
465  call check_status(nf90_get_var(ncid, variable_id, data1d))
466  else
467  ! 3D will reshape the data to take account of the column-row major C-F transposition
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)/))
471  deallocate(sdata)
472  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ read_specific_global_attribute()

character(len=:) function, allocatable, target kidreader_mod::read_specific_global_attribute ( integer, intent(in)  ncid,
character(len=*), intent(in)  key 
)
private

Will read a specific global NetCDF attribute.

Parameters
ncidThe NetCDF file id
keyThe name (key) of the attribute to read
Returns
String representation of the attribute value

Definition at line 521 of file kidreader.F90.

522  integer, intent(in) :: ncid
523  character(len=*), intent(in) :: key
524 
525  integer :: length
526  character(len=:),allocatable,target :: read_specific_global_attribute
527 
528  call check_status(nf90_inquire_attribute(ncid, nf90_global, key, len = length))
529  allocate(character(length) :: read_specific_global_attribute)
530  call check_status(nf90_get_att(ncid, nf90_global, key, read_specific_global_attribute))
Here is the call graph for this function:
Here is the caller graph for this function:

◆ read_variables()

subroutine kidreader_mod::read_variables ( integer, intent(in)  ncid,
integer, intent(in)  time_dim,
integer, intent(in)  z_dim,
integer, intent(in)  z_half_dim,
integer, intent(in)  x_dim,
integer, intent(in)  x_half_dim,
real, dimension(:), intent(inout), allocatable  time,
real, dimension(:), intent(inout), allocatable  x,
real, dimension(:), intent(inout), allocatable  z,
real, dimension(:), intent(inout), allocatable  x_half,
real, dimension(:), intent(inout), allocatable  z_half,
real, dimension(:,:,:), intent(inout), allocatable  u,
real, dimension(:,:,:), intent(inout), allocatable  w,
real, dimension(:,:,:), intent(inout), allocatable  v 
)
private

Reads the variables from the NetCDF KiD model file.

Parameters
ncidThe id of the NetCDF file
time_dimThe number of elements in the time dimension
z_dimThe number of elements in the z primal grid dimension
z_half_dimThe number of elements in the z dual grid dimension
x_dimThe number of elements in the x primal grid dimension
x_half_dimThe number of elements in the x dual grid dimension
timeThe time data field that is to be read
xThe points on the primal grid in x dimension
zThe points on the primal grid z dimension
x_halfThe points on the dual grid x dimension
z_halfThe points on the dual grid z dimension
uVelocity field in dimension x to read
wVelocity field in dimension z to read

Definition at line 404 of file kidreader.F90.

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
408 
409  allocate(time(time_dim))
410  allocate(z(z_dim))
411  allocate(x(x_dim))
412  allocate(z_half(z_half_dim))
413  allocate(x_half(x_half_dim))
414  ! Due to NetCDF being in C, need to reverse the data order for F2003
415  if (rotate_xy) then
416  allocate(v(time_dim, z_dim, x_half_dim))
417  else
418  allocate(u(time_dim, z_dim, x_half_dim))
419  end if
420  allocate(w(time_dim, z_half_dim, x_half_dim))
421 
422  call read_single_variable(ncid, time_key, data1d=time)
423  call read_single_variable(ncid, z_key, data1d=z)
424  call read_single_variable(ncid, z_half_key, data1d=z_half)
425  call read_single_variable(ncid, x_key, data1d=x)
426  call read_single_variable(ncid, x_half_key, data1d=x_half)
427  if (rotate_xy) then
428  call read_single_variable(ncid, u_key, data3d=v)
429  else
430  call read_single_variable(ncid, u_key, data3d=u)
431  end if
432  call read_single_variable(ncid, w_key, data3d=w)
433 
434  if (clone_to_3d) then
435  if (.not. allocated(v)) then
436  allocate(v(time_dim, z_dim, x_half_dim))
437  v=u
438  end if
439  if (.not. allocated(u)) then
440  allocate(u(time_dim, z_dim, x_half_dim))
441  u=v
442  end if
443  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ set_up_q_fields()

subroutine kidreader_mod::set_up_q_fields ( type(model_state_type), intent(inout)  current_state)
private

Definition at line 106 of file kidreader.F90.

107  type(model_state_type), intent(inout) :: current_state
108 
109  integer :: x_size, y_size, z_size, i
110 
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
114 
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
120  call initialise_single_q_field(current_state, i, z_size, y_size, x_size)
121  if (flood_q .and. i == 1) current_state%q(i)%data(:,:,:) = 1.
122  if (float_q .and. (.not. flood_q .or. i==2)) then
123  call populate_q_tracer(current_state, current_state%q(i))
124  end if
125  end do
Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ clone_to_3d

logical kidreader_mod::clone_to_3d = .false.
private

Definition at line 31 of file kidreader.F90.

◆ configuration_file

character(len=string_length) kidreader_mod::configuration_file
private

NetCDF model file to load.

Definition at line 37 of file kidreader.F90.

37  character(len=STRING_LENGTH) :: configuration_file

◆ domain_multiplication

integer kidreader_mod::domain_multiplication =1
private

The Q field coordinates configured by the user.

Definition at line 32 of file kidreader.F90.

32  integer :: domain_multiplication=1

◆ float_q

logical kidreader_mod::float_q = .false.
private

Definition at line 31 of file kidreader.F90.

◆ flood_q

logical kidreader_mod::flood_q = .false.
private

Definition at line 31 of file kidreader.F90.

31  logical :: flood_q = .false., float_q= .false., clone_to_3d = .false., rotate_xy=.false.

◆ number_q_coords

integer, parameter kidreader_mod::number_q_coords = 100
private

Number of Q field value coords that can be specified.

Definition at line 22 of file kidreader.F90.

22  integer, parameter :: NUMBER_Q_COORDS = 100

◆ q_coordinates_value

integer, dimension(number_q_coords), save kidreader_mod::q_coordinates_value
private

Definition at line 35 of file kidreader.F90.

◆ q_coordinates_x

integer, dimension(number_q_coords), save kidreader_mod::q_coordinates_x
private

Definition at line 35 of file kidreader.F90.

35  integer, dimension(NUMBER_Q_COORDS), save :: q_coordinates_x, q_coordinates_y, q_coordinates_z, q_coordinates_value

◆ q_coordinates_y

integer, dimension(number_q_coords), save kidreader_mod::q_coordinates_y
private

Definition at line 35 of file kidreader.F90.

◆ q_coordinates_z

integer, dimension(number_q_coords), save kidreader_mod::q_coordinates_z
private

Definition at line 35 of file kidreader.F90.

◆ rotate_xy

logical kidreader_mod::rotate_xy =.false.
private

Definition at line 31 of file kidreader.F90.

◆ time_key

character(len=*), parameter kidreader_mod::time_key = "time"
private

Corresponding NetCDF data time key.

Definition at line 23 of file kidreader.F90.

23  character(len=*), parameter :: TIME_KEY = "time",& !< Corresponding NetCDF data time key
24  z_key = "z",&
25  z_half_key = "z_half",&
26  x_key = "x",&
27  x_half_key = "x_half",&
28  u_key="u",&
29  w_key="w"

◆ u_key

character(len=*), parameter kidreader_mod::u_key ="u"
private

Corresponding NetCDF data u flow field key.

Definition at line 23 of file kidreader.F90.

◆ w_key

character(len=*), parameter kidreader_mod::w_key ="w"
private

Corresponding NetCDF data w flow field key.

Definition at line 23 of file kidreader.F90.

◆ x_half_key

character(len=*), parameter kidreader_mod::x_half_key = "x_half"
private

Corresponding NetCDF data x half (primal grid) key.

Definition at line 23 of file kidreader.F90.

◆ x_key

character(len=*), parameter kidreader_mod::x_key = "x"
private

Corresponding NetCDF data x (primal grid) key.

Definition at line 23 of file kidreader.F90.

◆ z_half_key

character(len=*), parameter kidreader_mod::z_half_key = "z_half"
private

Corresponding NetCDF data z half (dual grid) key.

Definition at line 23 of file kidreader.F90.

◆ z_key

character(len=*), parameter kidreader_mod::z_key = "z"
private

Corresponding NetCDF data z (primal grid) key.

Definition at line 23 of file kidreader.F90.

logging_mod::log_error
integer, parameter, public log_error
Only log ERROR messages.
Definition: logging.F90:11
logging_mod::log_info
integer, parameter, public log_info
Log INFO, WARNING and ERROR messages.
Definition: logging.F90:13
logging_mod::log_log
subroutine, public log_log(level, message, str)
Logs a message at the specified level. If the level is above the current level then the message is ig...
Definition: logging.F90:75
logging_mod::log_get_logging_level
integer function, public log_get_logging_level()
Retrieves the current logging level.
Definition: logging.F90:122
logging_mod::log_debug
integer, parameter, public log_debug
Log DEBUG, INFO, WARNING and ERROR messages.
Definition: logging.F90:14
logging_mod::log_master_log
subroutine, public log_master_log(level, message)
Will log just from the master process.
Definition: logging.F90:47
datadefn_mod::default_precision
integer, parameter, public default_precision
MPI communication type which we use for the prognostic and calculation data.
Definition: datadefn.F90:17