MONC
readcheckpoint.F90
Go to the documentation of this file.
1 
3 #ifndef TEST_MODE
4  use netcdf, only : nf90_global, nf90_nowrite, nf90_inquire_attribute, nf90_open, nf90_inq_dimid, nf90_inquire_dimension, &
5  nf90_inq_varid, nf90_get_var, nf90_get_att, nf90_close
6 #else
9 #endif
10  use datadefn_mod, only : string_length
11  use state_mod, only : model_state_type
26  implicit none
27 
28 #ifndef TEST_MODE
29  private
30 #endif
31 
33 
34 contains
35 
39  subroutine read_checkpoint_file(current_state, filename)
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.
71  end subroutine read_checkpoint_file
72 
76  subroutine decompose_grid(current_state)
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
84  end subroutine decompose_grid
85 
88  subroutine initalise_source_and_sav_fields(current_state)
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
130  end subroutine initalise_source_and_sav_fields
131 
135  subroutine load_misc(current_state, ncid)
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)
161  end subroutine load_misc
162 
167  function read_specific_global_attribute(ncid, key)
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))
177  end function read_specific_global_attribute
178 
188  subroutine load_all_fields(current_state, ncid)
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
257  end subroutine load_all_fields
258 
262  logical function does_field_exist(ncid, variable_key)
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)
269  end function does_field_exist
270 
273  subroutine load_q_indices(ncid)
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 /)))
287  call set_q_index(conv_to_integer(trim(value)), key)
288  end do
289  end if
290  end subroutine load_q_indices
291 
294  integer function get_number_q_indices(ncid)
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
306  end if
307  end function get_number_q_indices
308 
318  subroutine load_single_3d_field(ncid, local_grid, field, z_grid, y_grid, x_grid, variable_key, multi_process, fourth_dim_loc)
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.
373  end subroutine load_single_3d_field
374 
384  subroutine load_global_grid(current_state, ncid, z_dim, y_dim, x_dim, z_found, y_found, x_found)
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)
398  end subroutine load_global_grid
399 
400  subroutine load_mean_profiles(current_state, ncid, z_dim_id)
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
473  end subroutine load_mean_profiles
474 
480  subroutine define_vertical_levels(ncid, current_state, z_key, z_dim_id)
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)
503  end subroutine define_vertical_levels
504 
511  subroutine read_dimension_of_grid(ncid, grid, variable_key, dimension, dimension_id)
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.
537  end subroutine read_dimension_of_grid
538 
546  subroutine read_single_variable(ncid, key, int_data, real_data, real_data_1d, real_data_1d_double, real_data_2d_double, &
547  real_data_3d, integer_data_1d, start, count, map)
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
618  end subroutine read_single_variable
619 
632  subroutine read_dimensions(ncid, z_dim, y_dim, x_dim, z_found, y_found, x_found)
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)
640  end subroutine read_dimensions
character(len= *), parameter olubar
Generic add interface for adding different types of data to the databases.
character(len= *), parameter p_key
Pressure variable NetCDF key.
subroutine decompose_grid(current_state)
Calls out to the parallel decomposition strategy to decompose the grid based upon the configuration r...
character(len= *), parameter zn_key
character(len= *), parameter olzthbar
character(len= *), parameter y_key
character(len= *), parameter u_key
U variable NetCDF key.
character(len= *), parameter thref
integer function nf90_inquire_attribute(ncid, attributeid, key, len)
integer function nf90_get_att(ncid, attributeid, key, value)
subroutine, public read_checkpoint_file(current_state, filename)
Reads in a NetCDF checkpoint file and uses this to initialise the model.
type(standard_q_names_type), public standard_q_names
Definition: q_indices.F90:59
integer, parameter, public dual_grid
Definition: grids.F90:16
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.
character(len= *), parameter q_key
Q variable NetCDF key.
character(len= *), parameter empty_dim_key
Empty dimension key.
subroutine load_mean_profiles(current_state, ncid, z_dim_id)
character(len= *), parameter y_dim_key
Y dimension/variable key.
integer, parameter, public log_error
Only log ERROR messages.
Definition: logging.F90:11
Contains prognostic field definitions and functions.
Definition: prognostics.F90:2
character(len= *), parameter zv_key
A prognostic field which is assumed to be 3D.
Definition: prognostics.F90:13
integer function nf90_inq_dimid(ncid, key, dim_id)
Logging utility.
Definition: logging.F90:2
character(len= *), parameter olzqbar_anonymous_name
Common checkpoint functionality which is used by reader and writers to NetCDF checkpoints.
character(len= *), parameter created_attribute_key
character(len= *), parameter zq_field_anonymous_name
integer, parameter, public default_precision
MPI communication type which we use for the prognostic and calculation data.
Definition: datadefn.F90:17
character(len= *), parameter zth_key
integer, parameter, public z_index
Grid index parameters.
Definition: grids.F90:14
integer, parameter nf90_global
character(len= *), parameter x_key
integer, parameter max_string_length
Maximum string length (stored size)
character(len= *), parameter timestep
Timestep NetCDF key.
character(len= *), parameter time_key
character(len= *), parameter zu_key
character(len= *), parameter z_dim_key
Z dimension/variable key.
Contains common definitions for the data and datatypes used by MONC.
Definition: datadefn.F90:2
integer function get_number_q_indices(ncid)
Retrieve sthe number of Q indice entries in the NetCDF. This is optional, so may return 0...
character(len= *), parameter olzvbar
character(len= *), parameter v_key
V variable NetCDF key.
The ModelState which represents the current state of a run.
Definition: state.F90:39
subroutine, public log_master_log(level, message)
Will log just from the master process.
Definition: logging.F90:47
integer, parameter nf90_nowrite
Conversion between common inbuilt FORTRAN data types.
Definition: conversions.F90:5
integer function nf90_inq_varid(ncid, key, varid)
Converts data types to strings.
Definition: conversions.F90:36
character(len= *), parameter olvbar
character(len= *), parameter absolute_new_dtm_key
character(len= *), parameter q_indices_key
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 ha...
subroutine initalise_source_and_sav_fields(current_state)
Initialises the source and sav (for u,v,w) fields for allocated prognostics.
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
Converts data types to logical.
Definition: conversions.F90:69
type(q_metadata_type) function, public get_indices_descriptor(i)
Retrieves the indicies descriptor at a specific location.
Definition: q_indices.F90:100
character(len= *), parameter z_key
subroutine load_q_indices(ncid)
Loads the Q indices (if they exist) into the model Q index register.
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 ea...
Defines the global grid.
Definition: grids.F90:100
This manages the Q variables and specifically the mapping between names and the index that they are s...
Definition: q_indices.F90:2
character(len= *), parameter string_dim_key
String dimension key.
Determines whether a data item can be represented as a logical or not.
Definition: conversions.F90:98
character(len= *), parameter w_key
W variable NetCDF key.
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...
character(len= *), parameter title_attribute_key
character(len= *), parameter dtm_key
subroutine check_status(status, found_flag)
Will check a NetCDF status and write to log_log error any decoded statuses. Can be used to decode whe...
Defined the local grid, i.e. the grid held on this process after decomposition.
Definition: grids.F90:111
logical function does_field_exist(ncid, variable_key)
Determines whether a variable (field) exists within the NetCDF checkpoint file The NetCDF file id T...
character(len= *), parameter q_dim_key
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 ...
character(len= *), parameter olzubar
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 se...
integer function nf90_close(ncid)
Determines whether a data item can be represented as an integer or not.
Definition: conversions.F90:79
subroutine load_all_fields(current_state, ncid)
Reads in and initialises all prognostic fields. It will check the NetCDF checkpoint file and dependin...
character(len= *), parameter dtm_new_key
Converts data types to real.
Definition: conversions.F90:58
subroutine load_misc(current_state, ncid)
Loads in misc data from the checkpoint file.
integer, parameter, public string_length
Default length of strings.
Definition: datadefn.F90:10
Will read in a NetCDF checkpoint file and initialise the model state_mod based upon this...
character(len= *), parameter olzqbar
integer, parameter, public primal_grid
Grid type parameters (usually applied to each dimension of a prognostic)
Definition: grids.F90:16
character(len= *), parameter q_field_anonymous_name
character(len= *), parameter olthbar
Functionality to support the different types of grid and abstraction between global grids and local o...
Definition: grids.F90:5
Determines whether a data item can be represented as a real or not.
Definition: conversions.F90:89
character(len= *), parameter olqbar_anonymous_name
Manages the options database. Contains administration functions and deduce runtime options from the c...
integer, parameter, public log_info
Log INFO, WARNING and ERROR messages.
Definition: logging.F90:13
character(len= *), parameter q_indices_dim_key
character(len= *), parameter x_dim_key
X dimension/variable key.
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...
Converts data types to integers.
Definition: conversions.F90:47
character(len= *), parameter th_key
Theta variable NetCDF key.
character(len= *), parameter olqbar
character(len= *), parameter zq_key
character(len= *), parameter vgal
character(len= *), parameter zw_key
character(len= *), parameter ugal
The model state which represents the current state of a run.
Definition: state.F90:2
integer, parameter, public y_index
Definition: grids.F90:14
subroutine, public set_q_index(index, name)
Sets a Q index to be active at a specific index and sets the name.
Definition: q_indices.F90:71
character(len= *), parameter nqfields
integer, parameter, public x_index
Definition: grids.F90:14
integer function, public get_q_index(name, assigning_component)
Add in a new entry into the register if the name does not already exist or return the index of the pr...
Definition: q_indices.F90:112
integer function nf90_open(path, mode, ncid)
subroutine remove_null_terminator_from_string(net_cdf_string)
Removes NetCDF C style null termination of string. This is placed right at the end, after any spaces so trim will not actually trim any spaces due to null terminator.
integer function nf90_inquire_dimension(ncid, id, len)