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 @ncid The NetCDF file id @variable_key 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 load_pdf_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 78 of file readcheckpoint.F90.

79  type(model_state_type), intent(inout) :: current_state
80 
81  if (associated(current_state%parallel%decomposition_procedure)) then
82  call current_state%parallel%decomposition_procedure(current_state)
83  else
84  call log_log(log_error, "No decomposition specified")
85  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 511 of file readcheckpoint.F90.

512  type(model_state_type), intent(inout) :: current_state
513  character(len=*), intent(in) :: z_key
514  integer, intent(in) :: z_dim_id, ncid
515 
516  integer :: i, z_size
517  real, dimension(:), allocatable :: data
518 
519  call check_status(nf90_inquire_dimension(ncid, z_dim_id, len=z_size))
520  allocate(data(z_size))
521  call read_single_variable(ncid, z_key, real_data_1d=data)
522 
523  allocate(current_state%global_grid%configuration%vertical%kgd(z_size), &
524  current_state%global_grid%configuration%vertical%hgd(z_size), &
525  current_state%global_grid%configuration%vertical%thref(z_size))
526 
527  call read_single_variable(ncid, thref, real_data_1d_double=current_state%global_grid%configuration%vertical%thref)
528 
529  do i=1,z_size
530  current_state%global_grid%configuration%vertical%kgd(i) = i
531  current_state%global_grid%configuration%vertical%hgd(i) = real(data(i))
532  end do
533  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 @ncid The NetCDF file id @variable_key The NetCDF variable key we are looking up.

Definition at line 274 of file readcheckpoint.F90.

275  integer, intent(in) :: ncid
276  character(len=*), intent(in) :: variable_key
277 
278  integer :: variable_id
279 
280  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 306 of file readcheckpoint.F90.

307  integer, intent(in) :: ncid
308 
309  integer :: q_indices_dimid, q_indices_dim
310  logical :: found_flag
311 
312  call check_status(nf90_inq_dimid(ncid, q_indices_dim_key, q_indices_dimid), found_flag)
313  if (found_flag) then
314  call check_status(nf90_inquire_dimension(ncid, q_indices_dimid, len=q_indices_dim))
315  get_number_q_indices=q_indices_dim
316  else
317  get_number_q_indices=0
318  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 90 of file readcheckpoint.F90.

91  type(model_state_type), intent(inout) :: current_state
92 
93  integer :: x_size, y_size, z_size, i
94 
95  z_size = current_state%local_grid%size(z_index) + current_state%local_grid%halo_size(z_index) * 2
96  y_size = current_state%local_grid%size(y_index) + current_state%local_grid%halo_size(y_index) * 2
97  x_size = current_state%local_grid%size(x_index) + current_state%local_grid%halo_size(x_index) * 2
98 #ifdef U_ACTIVE
99  allocate(current_state%su%data(z_size, y_size, x_size))
100  current_state%su%data(:,:,:) = 0.
101  current_state%su%active=.true.
102  allocate(current_state%savu%data(z_size, y_size, x_size))
103  current_state%savu%data(:,:,:) = 0.
104  current_state%savu%active=.true.
105 #endif
106 #ifdef V_ACTIVE
107  allocate(current_state%sv%data(z_size, y_size, x_size))
108  current_state%sv%data(:,:,:) = 0.
109  current_state%sv%active=.true.
110  allocate(current_state%savv%data(z_size, y_size, x_size))
111  current_state%savv%data(:,:,:) = 0.
112  current_state%savv%active=.true.
113 #endif
114 #ifdef W_ACTIVE
115  allocate(current_state%sw%data(z_size, y_size, x_size))
116  current_state%sw%data(:,:,:) = 0.
117  current_state%sw%active=.true.
118  allocate(current_state%savw%data(z_size, y_size, x_size))
119  current_state%savw%data(:,:,:) = 0.
120  current_state%savw%active=.true.
121 #endif
122  if (current_state%th%active) then
123  allocate(current_state%sth%data(z_size, y_size, x_size))
124  current_state%sth%data(:,:,:) = 0.
125  current_state%sth%active=.true.
126  end if
127  do i=1,current_state%number_q_fields
128  current_state%sq(i)%active=.true.
129  allocate(current_state%sq(i)%data(z_size, y_size, x_size))
130  current_state%sq(i)%data(:,:,:) = 0.
131  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 192 of file readcheckpoint.F90.

193  type(model_state_type), intent(inout) :: current_state
194  integer, intent(in) :: ncid
195 
196  integer :: i
197  logical :: multi_process
198  type(q_metadata_type) :: q_metadata
199  character(len=STRING_LENGTH) :: q_field_name, zq_field_name
200 
201  multi_process = current_state%parallel%processes .gt. 1
202 
203  if (does_field_exist(ncid, u_key)) then
204  call load_single_3d_field(ncid, current_state%local_grid, current_state%u, dual_grid, &
205  dual_grid, primal_grid, u_key, multi_process)
206  call load_single_3d_field(ncid, current_state%local_grid, current_state%zu, dual_grid, &
207  dual_grid, primal_grid, zu_key, multi_process)
208  end if
209  if (does_field_exist(ncid, v_key)) then
210  call load_single_3d_field(ncid, current_state%local_grid, current_state%v, dual_grid, &
211  primal_grid, dual_grid, v_key, multi_process)
212  call load_single_3d_field(ncid, current_state%local_grid, current_state%zv, dual_grid, &
213  primal_grid, dual_grid, zv_key, multi_process)
214  end if
215  if (does_field_exist(ncid, w_key)) then
216  call load_single_3d_field(ncid, current_state%local_grid, current_state%w, primal_grid, &
217  dual_grid, dual_grid, w_key, multi_process)
218  call load_single_3d_field(ncid, current_state%local_grid, current_state%zw, primal_grid, &
219  dual_grid, dual_grid, zw_key, multi_process)
220  end if
221  if (does_field_exist(ncid, th_key)) then
222  call load_single_3d_field(ncid, current_state%local_grid, current_state%th, dual_grid, &
223  dual_grid, dual_grid, th_key, multi_process)
224  call load_single_3d_field(ncid, current_state%local_grid, current_state%zth, dual_grid, &
225  dual_grid, dual_grid, zth_key, multi_process)
226  end if
227  if (does_field_exist(ncid, p_key)) then
228  call load_single_3d_field(ncid, current_state%local_grid, current_state%p, dual_grid, &
229  dual_grid, dual_grid, p_key, multi_process)
230  end if
231  allocate(current_state%q(current_state%number_q_fields), current_state%zq(current_state%number_q_fields), &
232  current_state%sq(current_state%number_q_fields))
233  if (does_field_exist(ncid, q_key)) then
234  do i=1,current_state%number_q_fields
235  call load_single_3d_field(ncid, current_state%local_grid, current_state%q(i), dual_grid, &
236  dual_grid, dual_grid, q_key, multi_process, i)
237  call load_single_3d_field(ncid, current_state%local_grid, current_state%zq(i), dual_grid, &
238  dual_grid, dual_grid, zq_key, multi_process, i)
239  end do
240  else
241  do i=1,current_state%number_q_fields
242  q_metadata=get_indices_descriptor(i)
243  q_field_name=trim(q_key)//"_"//trim(q_metadata%name)
244  zq_field_name=trim(zq_key)//"_"//trim(q_metadata%name)
245  if (.not. does_field_exist(ncid, q_field_name)) then
246  q_field_name=trim(q_field_anonymous_name)//"_"//trim(conv_to_string(i))
247  zq_field_name=trim(zq_field_anonymous_name)//"_"//trim(conv_to_string(i))
248  if (.not. does_field_exist(ncid, q_field_name)) then
249  call log_log(log_error, "No entry in checkpoint file for Q field "//trim(conv_to_string(i)))
250  end if
251  end if
252  if (.not. does_field_exist(ncid, zq_field_name)) then
253  call log_log(log_error, "Missmatch between q and zq field name in the checkpoint file")
254  end if
255  call load_single_3d_field(ncid, current_state%local_grid, current_state%q(i), dual_grid, &
256  dual_grid, dual_grid, q_field_name, multi_process)
257  call load_single_3d_field(ncid, current_state%local_grid, current_state%zq(i), dual_grid, &
258  dual_grid, dual_grid, zq_field_name, multi_process)
259  end do
260  end if
261  if (does_field_exist(ncid, sth_lw_key)) then
262  call load_single_3d_field(ncid, current_state%local_grid, current_state%sth_lw, dual_grid, &
263  dual_grid, dual_grid, sth_lw_key, multi_process)
264  end if
265  if (does_field_exist(ncid, sth_lw_key)) then
266  call load_single_3d_field(ncid, current_state%local_grid, current_state%sth_sw, dual_grid, &
267  dual_grid, dual_grid, sth_sw_key, multi_process)
268  endif
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 396 of file readcheckpoint.F90.

397  type(model_state_type), intent(inout) :: current_state
398  integer, intent(in) :: ncid, z_dim, y_dim, x_dim
399  logical, intent(in) :: z_found, y_found, x_found
400 
401  current_state%global_grid%dimensions = 0
402 
403  if (z_found) then
404  call read_dimension_of_grid(ncid, current_state%global_grid, z_key, z_index, z_dim)
405  call define_vertical_levels(ncid, current_state, z_key, z_dim)
406 
407  end if
408  if (y_found) call read_dimension_of_grid(ncid, current_state%global_grid, y_key, y_index, y_dim)
409  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 412 of file readcheckpoint.F90.

413  type(model_state_type), intent(inout) :: current_state
414  integer, intent(in) :: ncid, z_dim_id
415 
416  integer :: z_size, i
417  type(q_metadata_type) :: q_metadata
418  character(len=STRING_LENGTH) :: q_field_name, zq_field_name
419 
420  call check_status(nf90_inquire_dimension(ncid, z_dim_id, len=z_size))
421  if (does_field_exist(ncid, olubar)) then
422  allocate(current_state%global_grid%configuration%vertical%olubar(z_size))
423  call read_single_variable(ncid, olubar, real_data_1d_double=current_state%global_grid%configuration%vertical%olubar)
424  end if
425  if (does_field_exist(ncid, olzubar)) then
426  allocate(current_state%global_grid%configuration%vertical%olzubar(z_size))
427  call read_single_variable(ncid, olzubar, real_data_1d_double=current_state%global_grid%configuration%vertical%olzubar)
428  end if
429  if (does_field_exist(ncid, olvbar)) then
430  allocate(current_state%global_grid%configuration%vertical%olvbar(z_size))
431  call read_single_variable(ncid, olvbar, real_data_1d_double=current_state%global_grid%configuration%vertical%olvbar)
432  end if
433  if (does_field_exist(ncid, olzvbar)) then
434  allocate(current_state%global_grid%configuration%vertical%olzvbar(z_size))
435  call read_single_variable(ncid, olzvbar, real_data_1d_double=current_state%global_grid%configuration%vertical%olzvbar)
436  end if
437  if (does_field_exist(ncid, olthbar)) then
438  allocate(current_state%global_grid%configuration%vertical%olthbar(z_size))
439  call read_single_variable(ncid, olthbar, real_data_1d_double=current_state%global_grid%configuration%vertical%olthbar)
440  end if
441  if (does_field_exist(ncid, olzthbar)) then
442  allocate(current_state%global_grid%configuration%vertical%olzthbar(z_size))
443  call read_single_variable(ncid, olzthbar, real_data_1d_double=current_state%global_grid%configuration%vertical%olzthbar)
444  end if
445  if (does_field_exist(ncid, olqbar)) then
446  allocate(current_state%global_grid%configuration%vertical%olqbar(z_size, current_state%number_q_fields))
447  call read_single_variable(ncid, olqbar, real_data_2d_double=current_state%global_grid%configuration%vertical%olqbar)
448  else if (current_state%number_q_fields .gt. 0) then
449  do i=1,current_state%number_q_fields
450  q_metadata=get_indices_descriptor(i)
451  q_field_name=trim(olqbar)//"_"//trim(q_metadata%name)
452  if (.not. does_field_exist(ncid, q_field_name)) then
453  q_field_name=trim(olqbar_anonymous_name)//"_"//trim(conv_to_string(i))
454  if (.not. does_field_exist(ncid, q_field_name)) then
455  cycle
456  end if
457  end if
458  if (.not. allocated(current_state%global_grid%configuration%vertical%olqbar)) then
459  allocate(current_state%global_grid%configuration%vertical%olqbar(z_size, current_state%number_q_fields))
460  end if
461  call read_single_variable(ncid, q_field_name, &
462  real_data_1d_double=current_state%global_grid%configuration%vertical%olqbar(:, i))
463  end do
464  end if
465  if (does_field_exist(ncid, olzqbar)) then
466  allocate(current_state%global_grid%configuration%vertical%olzqbar(z_size, current_state%number_q_fields))
467  call read_single_variable(ncid, olzqbar, real_data_2d_double=current_state%global_grid%configuration%vertical%olzqbar)
468  else if (current_state%number_q_fields .gt. 0) then
469  do i=1,current_state%number_q_fields
470  q_metadata=get_indices_descriptor(i)
471  q_field_name=trim(olzqbar)//"_"//trim(q_metadata%name)
472  if (.not. does_field_exist(ncid, q_field_name)) then
473  q_field_name=trim(olzqbar_anonymous_name)//"_"//trim(conv_to_string(i))
474  if (.not. does_field_exist(ncid, q_field_name)) then
475  cycle
476  end if
477  end if
478  if (.not. allocated(current_state%global_grid%configuration%vertical%olzqbar)) then
479  allocate(current_state%global_grid%configuration%vertical%olzqbar(z_size, current_state%number_q_fields))
480  end if
481  call read_single_variable(ncid, q_field_name, &
482  real_data_1d_double=current_state%global_grid%configuration%vertical%olzqbar(:, i))
483  end do
484  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 137 of file readcheckpoint.F90.

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

◆ load_pdf_profiles()

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

Definition at line 488 of file readcheckpoint.F90.

489  type(model_state_type), intent(inout) :: current_state
490  integer, intent(in) :: ncid, z_dim_id
491 
492  integer :: z_size, i
493 
494  call check_status(nf90_inquire_dimension(ncid, z_dim_id, len=z_size))
495  if (does_field_exist(ncid, wup)) then
496  allocate(current_state%global_grid%configuration%vertical%w_up(z_size))
497  call read_single_variable(ncid, wup, real_data_1d_double=current_state%global_grid%configuration%vertical%w_up)
498  end if
499  if (does_field_exist(ncid, wdwn)) then
500  allocate(current_state%global_grid%configuration%vertical%w_dwn(z_size))
501  call read_single_variable(ncid, wdwn, real_data_1d_double=current_state%global_grid%configuration%vertical%w_dwn)
502  end if
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 285 of file readcheckpoint.F90.

286  integer, intent(in) :: ncid
287 
288  integer :: number_q_indices, i, q_indices_id
289  character(len=MAX_STRING_LENGTH) :: key, value
290 
291  number_q_indices=get_number_q_indices(ncid)
292  if (number_q_indices .gt. 0) then
293  call check_status(nf90_inq_varid(ncid, q_indices_key, q_indices_id))
294  do i=1, number_q_indices
295  call check_status(nf90_get_var(ncid, q_indices_id, key, (/ 1, 1, i /)))
296  call check_status(nf90_get_var(ncid, q_indices_id, value, (/ 1, 2, i /)))
297  call remove_null_terminator_from_string(key)
298  call remove_null_terminator_from_string(value)
299  call set_q_index(conv_to_integer(trim(value)), key)
300  end do
301  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 330 of file readcheckpoint.F90.

331  character(len=*), intent(in) :: variable_key
332  integer, intent(in) :: z_grid, y_grid, x_grid, ncid
333  type(prognostic_field_type), intent(inout) :: field
334  type(local_grid_type), intent(inout) :: local_grid
335  logical, intent(in) :: multi_process
336  integer, optional, intent(in) :: fourth_dim_loc
337 
338  integer :: start(5), count(5), i, map(5)
339  integer :: variable_id, nd
340 
341  if (allocated(field%data)) deallocate(field%data)
342  allocate(field%data(local_grid%size(z_index) + local_grid%halo_size(z_index) * 2, local_grid%size(y_index) + &
343  local_grid%halo_size(y_index) * 2, local_grid%size(x_index) + local_grid%halo_size(x_index) * 2))
344  field%data(:,:,:)=0.
345 
346  if (multi_process .or. present(fourth_dim_loc)) then
347  do i=1,3
348  if (i==1) then
349  map(i)=1
350  else
351  map(i)=map(i-1)*local_grid%size(i-1)
352  end if
353  start(i) = local_grid%start(i)
354  count(i) = local_grid%size(i)
355  end do
356  if (present(fourth_dim_loc)) then
357  start(4) = fourth_dim_loc
358  count(4) = 1
359  map(4)=map(3)*local_grid%size(3)
360 
361  start(5)=1
362  count(5)=1
363  map(5)=map(4)
364  else
365  start(4)=1
366  map(4)=map(3)*local_grid%size(3)
367  count(4)=1
368  end if
369 
370  call read_single_variable(ncid, variable_key, real_data_3d=field%data(local_grid%local_domain_start_index(z_index):&
371  local_grid%local_domain_end_index(z_index),local_grid%local_domain_start_index(y_index):&
372  local_grid%local_domain_end_index(y_index), local_grid%local_domain_start_index(x_index):&
373  local_grid%local_domain_end_index(x_index)), start=start, count=count, map=map)
374  else
375  call read_single_variable(ncid, variable_key, real_data_3d=field%data(local_grid%local_domain_start_index(z_index):&
376  local_grid%local_domain_end_index(z_index),local_grid%local_domain_start_index(y_index):&
377  local_grid%local_domain_end_index(y_index), local_grid%local_domain_start_index(x_index):&
378  local_grid%local_domain_end_index(x_index)))
379  end if
380 
381  field%grid(z_index) = z_grid
382  field%grid(y_index) = y_grid
383  field%grid(x_index) = x_grid
384  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.

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

543  integer, intent(in) :: ncid, dimension_id, dimension
544  character(len=*), intent(in) :: variable_key
545  type(global_grid_type), intent(inout) :: grid
546 
547  integer :: dim_size
548 
549  call check_status(nf90_inquire_dimension(ncid, dimension_id, len=dim_size))
550 
551  if (variable_key .eq. "x" .or. variable_key .eq. "y") then
552  call read_single_variable(ncid, trim(variable_key)//"_resolution", real_data=grid%resolution(dimension))
553  call read_single_variable(ncid, trim(variable_key)//"_top", real_data=grid%top(dimension))
554  call read_single_variable(ncid, trim(variable_key)//"_bottom", real_data=grid%bottom(dimension))
555  else if (variable_key .eq. "z") then
556  allocate(grid%configuration%vertical%z(dim_size), grid%configuration%vertical%zn(dim_size))
557  call read_single_variable(ncid, z_key, real_data_1d_double=grid%configuration%vertical%z)
558  call read_single_variable(ncid, zn_key, real_data_1d_double=grid%configuration%vertical%zn)
559  grid%top(dimension) = int(grid%configuration%vertical%z(dim_size))
560  grid%resolution(dimension) = int(grid%configuration%vertical%z(2) - grid%configuration%vertical%z(1))
561  ! For now hard code the bottom of the grid in each dimension as 0, first element is the 0th + size
562  ! 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
563  grid%bottom(dimension) = 0
564  end if
565  grid%size(dimension) = dim_size
566  grid%dimensions = grid%dimensions + 1
567  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 663 of file readcheckpoint.F90.

664  integer, intent(in) :: ncid
665  integer, intent(out) :: z_dim, y_dim, x_dim
666  logical, intent(out) :: z_found, y_found, x_found
667 
668  call check_status(nf90_inq_dimid(ncid, z_dim_key, z_dim), z_found)
669  call check_status(nf90_inq_dimid(ncid, y_dim_key, y_dim), y_found)
670  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 577 of file readcheckpoint.F90.

579  integer, intent(in) :: ncid
580  character(len=*), intent(in) :: key
581  integer, intent(inout), optional :: int_data
582  real(kind=default_precision) , intent(inout), optional :: real_data
583  real, dimension(:), intent(inout), optional :: real_data_1d
584  real(kind=default_precision), dimension(:), intent(inout), optional :: real_data_1d_double
585  real(kind=default_precision), dimension(:,:), intent(inout), optional :: real_data_2d_double
586  real(kind=default_precision), dimension(:,:,:), intent(inout), optional :: real_data_3d
587  integer, dimension(:), intent(inout), optional :: integer_data_1d
588  integer, dimension(:), intent(in), optional :: start, count, map
589 
590  integer :: variable_id
591 
592  call check_status(nf90_inq_varid(ncid, key, variable_id))
593 
594  if (.not. present(int_data) .and. .not. present(real_data) .and. .not. present(real_data_1d) .and. &
595  .not. present(real_data_1d_double) .and. .not. present(real_data_2d_double) .and. &
596  .not. present(real_data_3d) .and. .not. present(integer_data_1d)) return
597 
598  if (present(int_data)) then
599  call check_status(nf90_get_var(ncid, variable_id, int_data))
600  return
601  end if
602 
603  if (present(real_data)) then
604  call check_status(nf90_get_var(ncid, variable_id, real_data))
605  return
606  end if
607 
608  if (present(real_data_1d)) then
609  if (present(start) .and. present(count) .and. present(map)) then
610  call check_status(nf90_get_var(ncid, variable_id, real_data_1d, start=start, count=count, map=map))
611  else if (present(start) .and. present(count)) then
612  call check_status(nf90_get_var(ncid, variable_id, real_data_1d, start=start, count=count))
613  else
614  call check_status(nf90_get_var(ncid, variable_id, real_data_1d))
615  end if
616  else if (present(real_data_1d_double)) then
617  if (present(start) .and. present(count) .and. present(map)) then
618  call check_status(nf90_get_var(ncid, variable_id, real_data_1d_double, start=start, count=count, map=map))
619  else if (present(start) .and. present(count)) then
620  call check_status(nf90_get_var(ncid, variable_id, real_data_1d_double, start=start, count=count))
621  else
622  call check_status(nf90_get_var(ncid, variable_id, real_data_1d_double))
623  end if
624  else if (present(real_data_2d_double)) then
625  if (present(start) .and. present(count) .and. present(map)) then
626  call check_status(nf90_get_var(ncid, variable_id, real_data_2d_double, start=start, count=count, map=map))
627  else if (present(start) .and. present(count)) then
628  call check_status(nf90_get_var(ncid, variable_id, real_data_2d_double, start=start, count=count))
629  else
630  call check_status(nf90_get_var(ncid, variable_id, real_data_2d_double))
631  end if
632  else if (present(real_data_3d)) then
633  if (present(start) .and. present(count) .and. present(map)) then
634  call check_status(nf90_get_var(ncid, variable_id, real_data_3d, start=start, count=count, map=map))
635  else if (present(start) .and. present(count)) then
636  call check_status(nf90_get_var(ncid, variable_id, real_data_3d, start=start, count=count))
637  else
638  call check_status(nf90_get_var(ncid, variable_id, real_data_3d))
639  end if
640  else if (present(integer_data_1d)) then
641  if (present(start) .and. present(count) .and. present(map)) then
642  call check_status(nf90_get_var(ncid, variable_id, integer_data_1d, start, count, map))
643  else if (present(start) .and. present(count)) then
644  call check_status(nf90_get_var(ncid, variable_id, integer_data_1d, start, count))
645  else
646  call check_status(nf90_get_var(ncid, variable_id, integer_data_1d))
647  end if
648  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 171 of file readcheckpoint.F90.

172  integer, intent(in) :: ncid
173  character(len=*), intent(in) :: key
174 
175  integer :: length
176  character(len=:),allocatable,target :: read_specific_global_attribute
177 
178  call check_status(nf90_inquire_attribute(ncid, nf90_global, key, len = length))
179  allocate(character(length) :: read_specific_global_attribute)
180  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:
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_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