MONC
Functions/Subroutines
checkpointer_read_checkpoint_mod Module Reference

Will read in a NetCDF checkpoint file and initialise the model state_mod based upon this. More...

Functions/Subroutines

subroutine, public read_checkpoint_file (current_state, filename)
 Reads in a NetCDF checkpoint file and uses this to initialise the model. More...
 
subroutine decompose_grid (current_state)
 Calls out to the parallel decomposition strategy to decompose the grid based upon the configuration read and number of processes. More...
 
subroutine initalise_source_and_sav_fields (current_state)
 Initialises the source and sav (for u,v,w) fields for allocated prognostics. More...
 
subroutine load_misc (current_state, ncid)
 Loads in misc data from the checkpoint file. More...
 
character(len=:) function, allocatable, target read_specific_global_attribute (ncid, key)
 Will read a global attribute from the checkpoint file - note that it allocates string memory here so the caller should deallocate this to avoid memory leaks. More...
 
subroutine load_all_fields (current_state, ncid)
 Reads in and initialises all prognostic fields. It will check the NetCDF checkpoint file and depending upon what data is in the file then it will set the corresponding fields up which in essence determines whether it has 1, 2 or 3D fields. This also supports 1, 2 or 3D grids_mod which will set any missing grid dimensions to be 1. More...
 
logical function does_field_exist (ncid, variable_key)
 Determines whether a variable (field) exists within the NetCDF checkpoint file The NetCDF file id The NetCDF variable key we are looking up. More...
 
subroutine load_q_indices (ncid)
 Loads the Q indices (if they exist) into the model Q index register. More...
 
integer function get_number_q_indices (ncid)
 Retrieve sthe number of Q indice entries in the NetCDF. This is optional, so may return 0. More...
 
subroutine load_single_3d_field (ncid, local_grid, field, z_grid, y_grid, x_grid, variable_key, multi_process, fourth_dim_loc)
 Loads in and initialises a single 3D prognostic field. Even if the model is being run in 1 or 2D, the fields are still stored in 3D with the missing dimension sizes being set to one. More...
 
subroutine load_global_grid (current_state, ncid, z_dim, y_dim, x_dim, z_found, y_found, x_found)
 Reads in and initialises the primal grid from the NetCDF checkpoint and will handle 1, 2 or 3D runs. More...
 
subroutine load_mean_profiles (current_state, ncid, z_dim_id)
 
subroutine define_vertical_levels (ncid, current_state, z_key, z_dim_id)
 Defines the vertical levels of the grid. This is both the grid points and corresponding height for each point in metres. More...
 
subroutine read_dimension_of_grid (ncid, grid, variable_key, dimension, dimension_id)
 Reads a specific dimension of the grid from the NetCDF and sets this up in the state_mod. More...
 
subroutine read_single_variable (ncid, key, int_data, real_data, real_data_1d, real_data_1d_double, real_data_2d_double, real_data_3d, integer_data_1d, start, count, map)
 Reads in a single variable and sets the values of the data based upon this. Handles reading in and setting vectors of 1 or 3D reals and 1d integers. More...
 
subroutine read_dimensions (ncid, z_dim, y_dim, x_dim, z_found, y_found, x_found)
 Will read in the dimension sizes from the NetCDF file along with whether the z, y and x dimensions have been found. If any have not been found then this corresponds to running the model with 2 or 1D grids_mod. More...
 

Detailed Description

Will read in a NetCDF checkpoint file and initialise the model state_mod based upon this.

Function/Subroutine Documentation

◆ decompose_grid()

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

Calls out to the parallel decomposition strategy to decompose the grid based upon the configuration read and number of processes.

Parameters
current_stateThe current model state

Definition at line 77 of file readcheckpoint.F90.

77  type(model_state_type), intent(inout) :: current_state
78 
79  if (associated(current_state%parallel%decomposition_procedure)) then
80  call current_state%parallel%decomposition_procedure(current_state)
81  else
82  call log_log(log_error, "No decomposition specified")
83  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ define_vertical_levels()

subroutine checkpointer_read_checkpoint_mod::define_vertical_levels ( integer, intent(in)  ncid,
type(model_state_type), intent(inout)  current_state,
character(len=*), intent(in)  z_key,
integer, intent(in)  z_dim_id 
)
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
z_keyKey of the Z grid points in the NetCDF
z_sizeNumber of grid points

Definition at line 481 of file readcheckpoint.F90.

481  type(model_state_type), intent(inout) :: current_state
482  character(len=*), intent(in) :: z_key
483  integer, intent(in) :: z_dim_id, ncid
484 
485  integer :: i, z_size
486  real, dimension(:), allocatable :: data
487 
488  call check_status(nf90_inquire_dimension(ncid, z_dim_id, len=z_size))
489  allocate(data(z_size))
490  call read_single_variable(ncid, z_key, real_data_1d=data)
491 
492  allocate(current_state%global_grid%configuration%vertical%kgd(z_size), &
493  current_state%global_grid%configuration%vertical%hgd(z_size), &
494  current_state%global_grid%configuration%vertical%thref(z_size))
495 
496  call read_single_variable(ncid, thref, real_data_1d_double=current_state%global_grid%configuration%vertical%thref)
497 
498  do i=1,z_size
499  current_state%global_grid%configuration%vertical%kgd(i) = i
500  current_state%global_grid%configuration%vertical%hgd(i) = real(data(i))
501  end do
502  deallocate(data)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ does_field_exist()

logical function checkpointer_read_checkpoint_mod::does_field_exist ( integer, intent(in)  ncid,
character(len=*), intent(in)  variable_key 
)
private

Determines whether a variable (field) exists within the NetCDF checkpoint file The NetCDF file id The NetCDF variable key we are looking up.

Definition at line 263 of file readcheckpoint.F90.

263  integer, intent(in) :: ncid
264  character(len=*), intent(in) :: variable_key
265 
266  integer :: variable_id
267 
268  call check_status(nf90_inq_varid(ncid, variable_key, variable_id), does_field_exist)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_number_q_indices()

integer function checkpointer_read_checkpoint_mod::get_number_q_indices ( integer, intent(in)  ncid)
private

Retrieve sthe number of Q indice entries in the NetCDF. This is optional, so may return 0.

Parameters
ncidThe NetCDF file id

Definition at line 295 of file readcheckpoint.F90.

295  integer, intent(in) :: ncid
296 
297  integer :: q_indices_dimid, q_indices_dim
298  logical :: found_flag
299 
300  call check_status(nf90_inq_dimid(ncid, q_indices_dim_key, q_indices_dimid), found_flag)
301  if (found_flag) then
302  call check_status(nf90_inquire_dimension(ncid, q_indices_dimid, len=q_indices_dim))
303  get_number_q_indices=q_indices_dim
304  else
305  get_number_q_indices=0
306  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ initalise_source_and_sav_fields()

subroutine checkpointer_read_checkpoint_mod::initalise_source_and_sav_fields ( type(model_state_type), intent(inout)  current_state)
private

Initialises the source and sav (for u,v,w) fields for allocated prognostics.

Parameters
current_stateThe model current state

Definition at line 89 of file readcheckpoint.F90.

89  type(model_state_type), intent(inout) :: current_state
90 
91  integer :: x_size, y_size, z_size, i
92 
93  z_size = current_state%local_grid%size(z_index) + current_state%local_grid%halo_size(z_index) * 2
94  y_size = current_state%local_grid%size(y_index) + current_state%local_grid%halo_size(y_index) * 2
95  x_size = current_state%local_grid%size(x_index) + current_state%local_grid%halo_size(x_index) * 2
96 #ifdef U_ACTIVE
97  allocate(current_state%su%data(z_size, y_size, x_size))
98  current_state%su%data(:,:,:) = 0.
99  current_state%su%active=.true.
100  allocate(current_state%savu%data(z_size, y_size, x_size))
101  current_state%savu%data(:,:,:) = 0.
102  current_state%savu%active=.true.
103 #endif
104 #ifdef V_ACTIVE
105  allocate(current_state%sv%data(z_size, y_size, x_size))
106  current_state%sv%data(:,:,:) = 0.
107  current_state%sv%active=.true.
108  allocate(current_state%savv%data(z_size, y_size, x_size))
109  current_state%savv%data(:,:,:) = 0.
110  current_state%savv%active=.true.
111 #endif
112 #ifdef W_ACTIVE
113  allocate(current_state%sw%data(z_size, y_size, x_size))
114  current_state%sw%data(:,:,:) = 0.
115  current_state%sw%active=.true.
116  allocate(current_state%savw%data(z_size, y_size, x_size))
117  current_state%savw%data(:,:,:) = 0.
118  current_state%savw%active=.true.
119 #endif
120  if (current_state%th%active) then
121  allocate(current_state%sth%data(z_size, y_size, x_size))
122  current_state%sth%data(:,:,:) = 0.
123  current_state%sth%active=.true.
124  end if
125  do i=1,current_state%number_q_fields
126  current_state%sq(i)%active=.true.
127  allocate(current_state%sq(i)%data(z_size, y_size, x_size))
128  current_state%sq(i)%data(:,:,:) = 0.
129  end do
Here is the caller graph for this function:

◆ load_all_fields()

subroutine checkpointer_read_checkpoint_mod::load_all_fields ( type(model_state_type), intent(inout)  current_state,
integer, intent(in)  ncid 
)
private

Reads in and initialises all prognostic fields. It will check the NetCDF checkpoint file and depending upon what data is in the file then it will set the corresponding fields up which in essence determines whether it has 1, 2 or 3D fields. This also supports 1, 2 or 3D grids_mod which will set any missing grid dimensions to be 1.

Parameters
current_stateThe current model state_mod
ncidThe NetCDF file id
z_dimThe size of the z dimension
y_dimThe size of the y dimension
x_dimThe size of the x dimension

Definition at line 189 of file readcheckpoint.F90.

189  type(model_state_type), intent(inout) :: current_state
190  integer, intent(in) :: ncid
191 
192  integer :: i
193  logical :: multi_process
194  type(q_metadata_type) :: q_metadata
195  character(len=STRING_LENGTH) :: q_field_name, zq_field_name
196 
197  multi_process = current_state%parallel%processes .gt. 1
198 
199  if (does_field_exist(ncid, u_key)) then
200  call load_single_3d_field(ncid, current_state%local_grid, current_state%u, dual_grid, &
201  dual_grid, primal_grid, u_key, multi_process)
202  call load_single_3d_field(ncid, current_state%local_grid, current_state%zu, dual_grid, &
203  dual_grid, primal_grid, zu_key, multi_process)
204  end if
205  if (does_field_exist(ncid, v_key)) then
206  call load_single_3d_field(ncid, current_state%local_grid, current_state%v, dual_grid, &
207  primal_grid, dual_grid, v_key, multi_process)
208  call load_single_3d_field(ncid, current_state%local_grid, current_state%zv, dual_grid, &
209  primal_grid, dual_grid, zv_key, multi_process)
210  end if
211  if (does_field_exist(ncid, w_key)) then
212  call load_single_3d_field(ncid, current_state%local_grid, current_state%w, primal_grid, &
213  dual_grid, dual_grid, w_key, multi_process)
214  call load_single_3d_field(ncid, current_state%local_grid, current_state%zw, primal_grid, &
215  dual_grid, dual_grid, zw_key, multi_process)
216  end if
217  if (does_field_exist(ncid, th_key)) then
218  call load_single_3d_field(ncid, current_state%local_grid, current_state%th, dual_grid, &
219  dual_grid, dual_grid, th_key, multi_process)
220  call load_single_3d_field(ncid, current_state%local_grid, current_state%zth, dual_grid, &
221  dual_grid, dual_grid, zth_key, multi_process)
222  end if
223  if (does_field_exist(ncid, p_key)) then
224  call load_single_3d_field(ncid, current_state%local_grid, current_state%p, dual_grid, &
225  dual_grid, dual_grid, p_key, multi_process)
226  end if
227  allocate(current_state%q(current_state%number_q_fields), current_state%zq(current_state%number_q_fields), &
228  current_state%sq(current_state%number_q_fields))
229  if (does_field_exist(ncid, q_key)) then
230  do i=1,current_state%number_q_fields
231  call load_single_3d_field(ncid, current_state%local_grid, current_state%q(i), dual_grid, &
232  dual_grid, dual_grid, q_key, multi_process, i)
233  call load_single_3d_field(ncid, current_state%local_grid, current_state%zq(i), dual_grid, &
234  dual_grid, dual_grid, zq_key, multi_process, i)
235  end do
236  else
237  do i=1,current_state%number_q_fields
238  q_metadata=get_indices_descriptor(i)
239  q_field_name=trim(q_key)//"_"//trim(q_metadata%name)
240  zq_field_name=trim(zq_key)//"_"//trim(q_metadata%name)
241  if (.not. does_field_exist(ncid, q_field_name)) then
242  q_field_name=trim(q_field_anonymous_name)//"_"//trim(conv_to_string(i))
243  zq_field_name=trim(zq_field_anonymous_name)//"_"//trim(conv_to_string(i))
244  if (.not. does_field_exist(ncid, q_field_name)) then
245  call log_log(log_error, "No entry in checkpoint file for Q field "//trim(conv_to_string(i)))
246  end if
247  end if
248  if (.not. does_field_exist(ncid, zq_field_name)) then
249  call log_log(log_error, "Missmatch between q and zq field name in the checkpoint file")
250  end if
251  call load_single_3d_field(ncid, current_state%local_grid, current_state%q(i), dual_grid, &
252  dual_grid, dual_grid, q_field_name, multi_process)
253  call load_single_3d_field(ncid, current_state%local_grid, current_state%zq(i), dual_grid, &
254  dual_grid, dual_grid, zq_field_name, multi_process)
255  end do
256  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ load_global_grid()

subroutine checkpointer_read_checkpoint_mod::load_global_grid ( type(model_state_type), intent(inout)  current_state,
integer, intent(in)  ncid,
integer, intent(in)  z_dim,
integer, intent(in)  y_dim,
integer, intent(in)  x_dim,
logical, intent(in)  z_found,
logical, intent(in)  y_found,
logical, intent(in)  x_found 
)
private

Reads in and initialises the primal grid from the NetCDF checkpoint and will handle 1, 2 or 3D runs.

Parameters
current_stateThe current model state_mod
ncidThe NetCDF file id
z_dimSize in the z dimension
y_dimSize in the y dimension
x_simSize in the x dimension
z_foundWhether the z grid dimension exists in the checkpoint
y_foundWhether the y grid dimension exists in the checkpoint
x_foundWhether the x grid dimension exists in the checkpoint

Definition at line 385 of file readcheckpoint.F90.

385  type(model_state_type), intent(inout) :: current_state
386  integer, intent(in) :: ncid, z_dim, y_dim, x_dim
387  logical, intent(in) :: z_found, y_found, x_found
388 
389  current_state%global_grid%dimensions = 0
390 
391  if (z_found) then
392  call read_dimension_of_grid(ncid, current_state%global_grid, z_key, z_index, z_dim)
393  call define_vertical_levels(ncid, current_state, z_key, z_dim)
394 
395  end if
396  if (y_found) call read_dimension_of_grid(ncid, current_state%global_grid, y_key, y_index, y_dim)
397  if (x_found) call read_dimension_of_grid(ncid, current_state%global_grid, x_key, x_index, x_dim)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ load_mean_profiles()

subroutine checkpointer_read_checkpoint_mod::load_mean_profiles ( type(model_state_type), intent(inout)  current_state,
integer, intent(in)  ncid,
integer, intent(in)  z_dim_id 
)
private

Definition at line 401 of file readcheckpoint.F90.

401  type(model_state_type), intent(inout) :: current_state
402  integer, intent(in) :: ncid, z_dim_id
403 
404  integer :: z_size, i
405  type(q_metadata_type) :: q_metadata
406  character(len=STRING_LENGTH) :: q_field_name, zq_field_name
407 
408  call check_status(nf90_inquire_dimension(ncid, z_dim_id, len=z_size))
409  if (does_field_exist(ncid, olubar)) then
410  allocate(current_state%global_grid%configuration%vertical%olubar(z_size))
411  call read_single_variable(ncid, olubar, real_data_1d_double=current_state%global_grid%configuration%vertical%olubar)
412  end if
413  if (does_field_exist(ncid, olzubar)) then
414  allocate(current_state%global_grid%configuration%vertical%olzubar(z_size))
415  call read_single_variable(ncid, olzubar, real_data_1d_double=current_state%global_grid%configuration%vertical%olzubar)
416  end if
417  if (does_field_exist(ncid, olvbar)) then
418  allocate(current_state%global_grid%configuration%vertical%olvbar(z_size))
419  call read_single_variable(ncid, olvbar, real_data_1d_double=current_state%global_grid%configuration%vertical%olvbar)
420  end if
421  if (does_field_exist(ncid, olzvbar)) then
422  allocate(current_state%global_grid%configuration%vertical%olzvbar(z_size))
423  call read_single_variable(ncid, olzvbar, real_data_1d_double=current_state%global_grid%configuration%vertical%olzvbar)
424  end if
425  if (does_field_exist(ncid, olthbar)) then
426  allocate(current_state%global_grid%configuration%vertical%olthbar(z_size))
427  call read_single_variable(ncid, olthbar, real_data_1d_double=current_state%global_grid%configuration%vertical%olthbar)
428  end if
429  if (does_field_exist(ncid, olzthbar)) then
430  allocate(current_state%global_grid%configuration%vertical%olzthbar(z_size))
431  call read_single_variable(ncid, olzthbar, real_data_1d_double=current_state%global_grid%configuration%vertical%olzthbar)
432  end if
433  if (does_field_exist(ncid, olqbar)) then
434  allocate(current_state%global_grid%configuration%vertical%olqbar(z_size, current_state%number_q_fields))
435  call read_single_variable(ncid, olqbar, real_data_2d_double=current_state%global_grid%configuration%vertical%olqbar)
436  else if (current_state%number_q_fields .gt. 0) then
437  do i=1,current_state%number_q_fields
438  q_metadata=get_indices_descriptor(i)
439  q_field_name=trim(olqbar)//"_"//trim(q_metadata%name)
440  if (.not. does_field_exist(ncid, q_field_name)) then
441  q_field_name=trim(olqbar_anonymous_name)//"_"//trim(conv_to_string(i))
442  if (.not. does_field_exist(ncid, q_field_name)) then
443  cycle
444  end if
445  end if
446  if (.not. allocated(current_state%global_grid%configuration%vertical%olqbar)) then
447  allocate(current_state%global_grid%configuration%vertical%olqbar(z_size, current_state%number_q_fields))
448  end if
449  call read_single_variable(ncid, q_field_name, &
450  real_data_1d_double=current_state%global_grid%configuration%vertical%olqbar(:, i))
451  end do
452  end if
453  if (does_field_exist(ncid, olzqbar)) then
454  allocate(current_state%global_grid%configuration%vertical%olzqbar(z_size, current_state%number_q_fields))
455  call read_single_variable(ncid, olzqbar, real_data_2d_double=current_state%global_grid%configuration%vertical%olzqbar)
456  else if (current_state%number_q_fields .gt. 0) then
457  do i=1,current_state%number_q_fields
458  q_metadata=get_indices_descriptor(i)
459  q_field_name=trim(olzqbar)//"_"//trim(q_metadata%name)
460  if (.not. does_field_exist(ncid, q_field_name)) then
461  q_field_name=trim(olzqbar_anonymous_name)//"_"//trim(conv_to_string(i))
462  if (.not. does_field_exist(ncid, q_field_name)) then
463  cycle
464  end if
465  end if
466  if (.not. allocated(current_state%global_grid%configuration%vertical%olzqbar)) then
467  allocate(current_state%global_grid%configuration%vertical%olzqbar(z_size, current_state%number_q_fields))
468  end if
469  call read_single_variable(ncid, q_field_name, &
470  real_data_1d_double=current_state%global_grid%configuration%vertical%olzqbar(:, i))
471  end do
472  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ load_misc()

subroutine checkpointer_read_checkpoint_mod::load_misc ( type(model_state_type), intent(inout)  current_state,
integer, intent(in)  ncid 
)
private

Loads in misc data from the checkpoint file.

Parameters
current_stateThe current model state_mod
ncidThe NetCDF file id

Definition at line 136 of file readcheckpoint.F90.

136  type(model_state_type), intent(inout) :: current_state
137  integer, intent(in) :: ncid
138 
139  integer :: i_data(1) ! Procedure requires a vector rather than scalar
140  real(kind=DEFAULT_PRECISION) :: r_data(1)
141 
142  call read_single_variable(ncid, timestep, integer_data_1d=i_data)
143  current_state%timestep = i_data(1)+1 ! plus one to increment for next timestep
144  !current_state%start_timestep = current_state%timestep
145  call read_single_variable(ncid, ugal, real_data_1d_double=r_data)
146  current_state%ugal = r_data(1)
147  call read_single_variable(ncid, vgal, real_data_1d_double=r_data)
148  current_state%vgal = r_data(1)
149  call read_single_variable(ncid, nqfields, integer_data_1d=i_data)
150  current_state%number_q_fields = i_data(1)
151  call read_single_variable(ncid, dtm_key, real_data_1d_double=r_data)
152  current_state%dtm = r_data(1)
153  call read_single_variable(ncid, dtm_new_key, real_data_1d_double=r_data)
154  current_state%dtm_new = r_data(1)
155  current_state%update_dtm = current_state%dtm .ne. current_state%dtm_new
156  call read_single_variable(ncid, absolute_new_dtm_key, real_data_1d_double=r_data)
157  current_state%absolute_new_dtm = r_data(1)
158  call read_single_variable(ncid, time_key, real_data_1d_double=r_data)
159  ! The time is written into checkpoint as time+dtm, therefore the time as read in has been correctly advanced
160  current_state%time = r_data(1)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ load_q_indices()

subroutine checkpointer_read_checkpoint_mod::load_q_indices ( integer, intent(in)  ncid)
private

Loads the Q indices (if they exist) into the model Q index register.

Parameters
ncidThe NetCDF file id

Definition at line 274 of file readcheckpoint.F90.

274  integer, intent(in) :: ncid
275 
276  integer :: number_q_indices, i, q_indices_id
277  character(len=MAX_STRING_LENGTH) :: key, value
278 
279  number_q_indices=get_number_q_indices(ncid)
280  if (number_q_indices .gt. 0) then
281  call check_status(nf90_inq_varid(ncid, q_indices_key, q_indices_id))
282  do i=1, number_q_indices
283  call check_status(nf90_get_var(ncid, q_indices_id, key, (/ 1, 1, i /)))
284  call check_status(nf90_get_var(ncid, q_indices_id, value, (/ 1, 2, i /)))
285  call remove_null_terminator_from_string(key)
286  call remove_null_terminator_from_string(value)
287  call set_q_index(conv_to_integer(trim(value)), key)
288  end do
289  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ load_single_3d_field()

subroutine checkpointer_read_checkpoint_mod::load_single_3d_field ( integer, intent(in)  ncid,
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,
character(len=*), intent(in)  variable_key,
logical, intent(in)  multi_process,
integer, intent(in), optional  fourth_dim_loc 
)
private

Loads in and initialises a single 3D prognostic field. Even if the model is being run in 1 or 2D, the fields are still stored in 3D with the missing dimension sizes being set to one.

Parameters
ncidThe NetCDF file id
fieldThe prognostic field that will be initialised from the checkpoint
gridThe grid that the prognostic field is on
variable_keyThe NetCDF variable name
dim_oneSize in dimension one
dim_twoSize in dimension two
dim_threeSize in dimension three

Definition at line 319 of file readcheckpoint.F90.

319  character(len=*), intent(in) :: variable_key
320  integer, intent(in) :: z_grid, y_grid, x_grid, ncid
321  type(prognostic_field_type), intent(inout) :: field
322  type(local_grid_type), intent(inout) :: local_grid
323  logical, intent(in) :: multi_process
324  integer, optional, intent(in) :: fourth_dim_loc
325 
326  integer :: start(5), count(5), i, map(5)
327  integer :: variable_id, nd
328 
329  if (allocated(field%data)) deallocate(field%data)
330  allocate(field%data(local_grid%size(z_index) + local_grid%halo_size(z_index) * 2, local_grid%size(y_index) + &
331  local_grid%halo_size(y_index) * 2, local_grid%size(x_index) + local_grid%halo_size(x_index) * 2))
332  field%data(:,:,:)=0.
333 
334  if (multi_process .or. present(fourth_dim_loc)) then
335  do i=1,3
336  if (i==1) then
337  map(i)=1
338  else
339  map(i)=map(i-1)*local_grid%size(i-1)
340  end if
341  start(i) = local_grid%start(i)
342  count(i) = local_grid%size(i)
343  end do
344  if (present(fourth_dim_loc)) then
345  start(4) = fourth_dim_loc
346  count(4) = 1
347  map(4)=map(3)*local_grid%size(3)
348 
349  start(5)=1
350  count(5)=1
351  map(5)=map(4)
352  else
353  start(4)=1
354  map(4)=map(3)*local_grid%size(3)
355  count(4)=1
356  end if
357 
358  call read_single_variable(ncid, variable_key, real_data_3d=field%data(local_grid%local_domain_start_index(z_index):&
359  local_grid%local_domain_end_index(z_index),local_grid%local_domain_start_index(y_index):&
360  local_grid%local_domain_end_index(y_index), local_grid%local_domain_start_index(x_index):&
361  local_grid%local_domain_end_index(x_index)), start=start, count=count, map=map)
362  else
363  call read_single_variable(ncid, variable_key, real_data_3d=field%data(local_grid%local_domain_start_index(z_index):&
364  local_grid%local_domain_end_index(z_index),local_grid%local_domain_start_index(y_index):&
365  local_grid%local_domain_end_index(y_index), local_grid%local_domain_start_index(x_index):&
366  local_grid%local_domain_end_index(x_index)))
367  end if
368 
369  field%grid(z_index) = z_grid
370  field%grid(y_index) = y_grid
371  field%grid(x_index) = x_grid
372  field%active = .true.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ read_checkpoint_file()

subroutine, public checkpointer_read_checkpoint_mod::read_checkpoint_file ( type(model_state_type), intent(inout)  current_state,
character(len=*), intent(in)  filename 
)

Reads in a NetCDF checkpoint file and uses this to initialise the model.

Parameters
currentStateThe current model state_mod
filenameThe filename of the checkpoint file to load

Definition at line 40 of file readcheckpoint.F90.

40  type(model_state_type), intent(inout) :: current_state
41  character(len=*), intent(in) :: filename
42 
43  integer :: ncid, z_dim, y_dim, x_dim
44  logical :: z_found, y_found, x_found
45  character(len=:), allocatable :: attribute_value
46 
47  call check_status(nf90_open(path = filename, mode = nf90_nowrite, ncid = ncid))
48  attribute_value=read_specific_global_attribute(ncid, "created")
49  call read_dimensions(ncid, z_dim, y_dim, x_dim, z_found, y_found, x_found)
50  call load_global_grid(current_state, ncid, z_dim, y_dim, x_dim, z_found, y_found, x_found)
51 
52  call decompose_grid(current_state)
53 
54  call load_q_indices(ncid)
55  call load_misc(current_state, ncid)
56  call load_all_fields(current_state, ncid)
57  call initalise_source_and_sav_fields(current_state)
58  call load_mean_profiles(current_state, ncid, z_dim)
59  call check_status(nf90_close(ncid))
60 
61  if (current_state%number_q_fields .gt. 0) then
62  ! Retrieve standard q indices
63  current_state%water_vapour_mixing_ratio_index=get_q_index(standard_q_names%VAPOUR, 'checkpoint')
64  current_state%liquid_water_mixing_ratio_index=get_q_index(standard_q_names%CLOUD_LIQUID_MASS, 'checkpoint')
65  end if
66 
67  call log_master_log(log_info, "Restarted configuration from checkpoint file `"//trim(filename)//"` created at "&
68  //attribute_value)
69  deallocate(attribute_value)
70  current_state%initialised=.true.
Here is the call graph for this function:

◆ read_dimension_of_grid()

subroutine checkpointer_read_checkpoint_mod::read_dimension_of_grid ( integer, intent(in)  ncid,
type(global_grid_type), intent(inout)  grid,
character(len=*), intent(in)  variable_key,
integer, intent(in)  dimension,
integer, intent(in)  dimension_id 
)
private

Reads a specific dimension of the grid from the NetCDF and sets this up in the state_mod.

Parameters
ncidThe NetCDF file id
gridThe model grid that this dimension lies on
variable_keyThe NetCDF variable name that we are reading
dimensionThe dimension corresponding to the definition in the grids_mod module
dimension_sizeNumber of grid points in this dimension

Definition at line 512 of file readcheckpoint.F90.

512  integer, intent(in) :: ncid, dimension_id, dimension
513  character(len=*), intent(in) :: variable_key
514  type(global_grid_type), intent(inout) :: grid
515 
516  integer :: dim_size
517 
518  call check_status(nf90_inquire_dimension(ncid, dimension_id, len=dim_size))
519 
520  if (variable_key .eq. "x" .or. variable_key .eq. "y") then
521  call read_single_variable(ncid, trim(variable_key)//"_resolution", real_data=grid%resolution(dimension))
522  call read_single_variable(ncid, trim(variable_key)//"_top", real_data=grid%top(dimension))
523  call read_single_variable(ncid, trim(variable_key)//"_bottom", real_data=grid%bottom(dimension))
524  else if (variable_key .eq. "z") then
525  allocate(grid%configuration%vertical%z(dim_size), grid%configuration%vertical%zn(dim_size))
526  call read_single_variable(ncid, z_key, real_data_1d_double=grid%configuration%vertical%z)
527  call read_single_variable(ncid, zn_key, real_data_1d_double=grid%configuration%vertical%zn)
528  grid%top(dimension) = int(grid%configuration%vertical%z(dim_size))
529  grid%resolution(dimension) = int(grid%configuration%vertical%z(2) - grid%configuration%vertical%z(1))
530  ! For now hard code the bottom of the grid in each dimension as 0, first element is the 0th + size
531  ! I.e. if dxx=100, start=0 then first point is 100 rather than 0 (in x and y), is correct in z in CP file
532  grid%bottom(dimension) = 0
533  end if
534  grid%size(dimension) = dim_size
535  grid%dimensions = grid%dimensions + 1
536  grid%active(dimension) = .true.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ read_dimensions()

subroutine checkpointer_read_checkpoint_mod::read_dimensions ( integer, intent(in)  ncid,
integer, intent(out)  z_dim,
integer, intent(out)  y_dim,
integer, intent(out)  x_dim,
logical, intent(out)  z_found,
logical, intent(out)  y_found,
logical, intent(out)  x_found 
)
private

Will read in the dimension sizes from the NetCDF file along with whether the z, y and x dimensions have been found. If any have not been found then this corresponds to running the model with 2 or 1D grids_mod.

Parameters
ncidThe NetCDf file id
z_dimSize in the z dimension (set to size of empty dimension if not found)
y_dimSize in the y dimension (set to size of empty dimension if not found)
x_dimSize in the x dimension (set to size of empty dimension if not found)
string_dimSize of the string dimension (size that we use to allocate strings)
key_value_pair_dimSize of key-value pair dimension which should be two
options_dimSize of the options dimension which corresponds to how many options have been provided
z_foundWhether the z dimension exists (grid exists in z dimension)
y_foundWhether the y dimension exists (grid exists in y dimension)
x_foundWhether the x dimension exists (grid exists in x dimension)

Definition at line 633 of file readcheckpoint.F90.

633  integer, intent(in) :: ncid
634  integer, intent(out) :: z_dim, y_dim, x_dim
635  logical, intent(out) :: z_found, y_found, x_found
636 
637  call check_status(nf90_inq_dimid(ncid, z_dim_key, z_dim), z_found)
638  call check_status(nf90_inq_dimid(ncid, y_dim_key, y_dim), y_found)
639  call check_status(nf90_inq_dimid(ncid, x_dim_key, x_dim), x_found)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ read_single_variable()

subroutine checkpointer_read_checkpoint_mod::read_single_variable ( integer, intent(in)  ncid,
character(len=*), intent(in)  key,
integer, intent(inout), optional  int_data,
real(kind=default_precision), intent(inout), optional  real_data,
real, dimension(:), intent(inout), optional  real_data_1d,
real(kind=default_precision), dimension(:), intent(inout), optional  real_data_1d_double,
real(kind=default_precision), dimension(:,:), intent(inout), optional  real_data_2d_double,
real(kind=default_precision), dimension(:,:,:), intent(inout), optional  real_data_3d,
integer, dimension(:), intent(inout), optional  integer_data_1d,
integer, dimension(:), intent(in), optional  start,
integer, dimension(:), intent(in), optional  count,
integer, dimension(:), intent(in), optional  map 
)
private

Reads in a single variable and sets the values of the data based upon this. Handles reading in and setting vectors of 1 or 3D reals and 1d integers.

Parameters
ncidThe NetCDF file id
keyThe variable key to read in
realData1DVector of 1D real data to read (optional)
realData3DVector of 3D real data to read (optional)
integerData1DVector of 1D integer data to read (optional)

Definition at line 548 of file readcheckpoint.F90.

548  integer, intent(in) :: ncid
549  character(len=*), intent(in) :: key
550  integer, intent(inout), optional :: int_data
551  real(kind=DEFAULT_PRECISION) , intent(inout), optional :: real_data
552  real, dimension(:), intent(inout), optional :: real_data_1d
553  real(kind=DEFAULT_PRECISION), dimension(:), intent(inout), optional :: real_data_1d_double
554  real(kind=DEFAULT_PRECISION), dimension(:,:), intent(inout), optional :: real_data_2d_double
555  real(kind=DEFAULT_PRECISION), dimension(:,:,:), intent(inout), optional :: real_data_3d
556  integer, dimension(:), intent(inout), optional :: integer_data_1d
557  integer, dimension(:), intent(in), optional :: start, count, map
558 
559  integer :: variable_id
560 
561  call check_status(nf90_inq_varid(ncid, key, variable_id))
562 
563  if (.not. present(int_data) .and. .not. present(real_data) .and. .not. present(real_data_1d) .and. &
564  .not. present(real_data_1d_double) .and. .not. present(real_data_2d_double) .and. &
565  .not. present(real_data_3d) .and. .not. present(integer_data_1d)) return
566 
567  if (present(int_data)) then
568  call check_status(nf90_get_var(ncid, variable_id, int_data))
569  return
570  end if
571 
572  if (present(real_data)) then
573  call check_status(nf90_get_var(ncid, variable_id, real_data))
574  return
575  end if
576 
577  if (present(real_data_1d)) then
578  if (present(start) .and. present(count) .and. present(map)) then
579  call check_status(nf90_get_var(ncid, variable_id, real_data_1d, start=start, count=count, map=map))
580  else if (present(start) .and. present(count)) then
581  call check_status(nf90_get_var(ncid, variable_id, real_data_1d, start=start, count=count))
582  else
583  call check_status(nf90_get_var(ncid, variable_id, real_data_1d))
584  end if
585  else if (present(real_data_1d_double)) then
586  if (present(start) .and. present(count) .and. present(map)) then
587  call check_status(nf90_get_var(ncid, variable_id, real_data_1d_double, start=start, count=count, map=map))
588  else if (present(start) .and. present(count)) then
589  call check_status(nf90_get_var(ncid, variable_id, real_data_1d_double, start=start, count=count))
590  else
591  call check_status(nf90_get_var(ncid, variable_id, real_data_1d_double))
592  end if
593  else if (present(real_data_2d_double)) then
594  if (present(start) .and. present(count) .and. present(map)) then
595  call check_status(nf90_get_var(ncid, variable_id, real_data_2d_double, start=start, count=count, map=map))
596  else if (present(start) .and. present(count)) then
597  call check_status(nf90_get_var(ncid, variable_id, real_data_2d_double, start=start, count=count))
598  else
599  call check_status(nf90_get_var(ncid, variable_id, real_data_2d_double))
600  end if
601  else if (present(real_data_3d)) then
602  if (present(start) .and. present(count) .and. present(map)) then
603  call check_status(nf90_get_var(ncid, variable_id, real_data_3d, start=start, count=count, map=map))
604  else if (present(start) .and. present(count)) then
605  call check_status(nf90_get_var(ncid, variable_id, real_data_3d, start=start, count=count))
606  else
607  call check_status(nf90_get_var(ncid, variable_id, real_data_3d))
608  end if
609  else if (present(integer_data_1d)) then
610  if (present(start) .and. present(count) .and. present(map)) then
611  call check_status(nf90_get_var(ncid, variable_id, integer_data_1d, start, count, map))
612  else if (present(start) .and. present(count)) then
613  call check_status(nf90_get_var(ncid, variable_id, integer_data_1d, start, count))
614  else
615  call check_status(nf90_get_var(ncid, variable_id, integer_data_1d))
616  end if
617  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 checkpointer_read_checkpoint_mod::read_specific_global_attribute ( integer, intent(in)  ncid,
character(len=*), intent(in)  key 
)
private

Will read a global attribute from the checkpoint file - note that it allocates string memory here so the caller should deallocate this to avoid memory leaks.

Parameters
ncidThe NetCDF file id
keyThe NetCDF global attributes key to look up

Definition at line 168 of file readcheckpoint.F90.

168  integer, intent(in) :: ncid
169  character(len=*), intent(in) :: key
170 
171  integer :: length
172  character(len=:),allocatable,target :: read_specific_global_attribute
173 
174  call check_status(nf90_inquire_attribute(ncid, nf90_global, key, len = length))
175  allocate(character(length) :: read_specific_global_attribute)
176  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: