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
19 z_dim_key,
q_dim_key,
q_key,
zq_key,
th_key,
zth_key,
p_key,
u_key,
v_key,
w_key,
zu_key,
zv_key,
zw_key,
x_key,
y_key, &
20 z_key,
zn_key,
nqfields,
ugal,
vgal,
time_key,
timestep,
created_attribute_key,
title_attribute_key,
absolute_new_dtm_key, &
22 max_string_length,
thref,
olubar,
olzubar,
olvbar,
olzvbar,
olthbar,
olzthbar,
olqbar,
olzqbar,
olqbar_anonymous_name, &
41 character(len=*),
intent(in) :: filename
43 integer :: ncid, z_dim, y_dim, x_dim
44 logical :: z_found, y_found, x_found
45 character(len=:),
allocatable :: attribute_value
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)
61 if (current_state%number_q_fields .gt. 0)
then 67 call log_master_log(
log_info,
"Restarted configuration from checkpoint file `"//trim(filename)//
"` created at "&
69 deallocate(attribute_value)
70 current_state%initialised=.true.
79 if (
associated(current_state%parallel%decomposition_procedure))
then 80 call current_state%parallel%decomposition_procedure(current_state)
91 integer :: x_size, y_size, z_size, i
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
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.
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.
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.
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.
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.
135 subroutine load_misc(current_state, ncid)
137 integer,
intent(in) :: ncid
140 real(kind=DEFAULT_PRECISION) :: r_data(1)
143 current_state%timestep = i_data(1)+1
146 current_state%ugal = r_data(1)
148 current_state%vgal = r_data(1)
150 current_state%number_q_fields = i_data(1)
152 current_state%dtm = r_data(1)
154 current_state%dtm_new = r_data(1)
155 current_state%update_dtm = current_state%dtm .ne. current_state%dtm_new
157 current_state%absolute_new_dtm = r_data(1)
160 current_state%time = r_data(1)
168 integer,
intent(in) :: ncid
169 character(len=*),
intent(in) :: key
172 character(len=:),
allocatable,
target :: read_specific_global_attribute
175 allocate(
character(length) :: read_specific_global_attribute)
190 integer,
intent(in) :: ncid
193 logical :: multi_process
195 character(len=STRING_LENGTH) :: q_field_name, zq_field_name
197 multi_process = current_state%parallel%processes .gt. 1
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))
230 do i=1,current_state%number_q_fields
237 do i=1,current_state%number_q_fields
239 q_field_name=trim(
q_key)//
"_"//trim(q_metadata%name)
240 zq_field_name=trim(
zq_key)//
"_"//trim(q_metadata%name)
249 call log_log(
log_error,
"Missmatch between q and zq field name in the checkpoint file")
263 integer,
intent(in) :: ncid
264 character(len=*),
intent(in) :: variable_key
266 integer :: variable_id
274 integer,
intent(in) :: ncid
276 integer :: number_q_indices, i, q_indices_id
277 character(len=MAX_STRING_LENGTH) :: key, value
280 if (number_q_indices .gt. 0)
then 282 do i=1, number_q_indices
295 integer,
intent(in) :: ncid
297 integer :: q_indices_dimid, q_indices_dim
298 logical :: found_flag
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
323 logical,
intent(in) :: multi_process
324 integer,
optional,
intent(in) :: fourth_dim_loc
326 integer :: start(5), count(5), i, map(5)
327 integer :: variable_id, nd
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) + &
334 if (multi_process .or.
present(fourth_dim_loc))
then 339 map(i)=map(i-1)*local_grid%size(i-1)
341 start(i) = local_grid%start(i)
342 count(i) = local_grid%size(i)
344 if (
present(fourth_dim_loc))
then 345 start(4) = fourth_dim_loc
347 map(4)=map(3)*local_grid%size(3)
354 map(4)=map(3)*local_grid%size(3)
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)
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)))
372 field%active = .true.
384 subroutine load_global_grid(current_state, ncid, z_dim, y_dim, x_dim, z_found, y_found, x_found)
386 integer,
intent(in) :: ncid, z_dim, y_dim, x_dim
387 logical,
intent(in) :: z_found, y_found, x_found
389 current_state%global_grid%dimensions = 0
402 integer,
intent(in) :: ncid, z_dim_id
406 character(len=STRING_LENGTH) :: q_field_name, zq_field_name
410 allocate(current_state%global_grid%configuration%vertical%olubar(z_size))
414 allocate(current_state%global_grid%configuration%vertical%olzubar(z_size))
418 allocate(current_state%global_grid%configuration%vertical%olvbar(z_size))
422 allocate(current_state%global_grid%configuration%vertical%olzvbar(z_size))
426 allocate(current_state%global_grid%configuration%vertical%olthbar(z_size))
430 allocate(current_state%global_grid%configuration%vertical%olzthbar(z_size))
434 allocate(current_state%global_grid%configuration%vertical%olqbar(z_size, current_state%number_q_fields))
436 else if (current_state%number_q_fields .gt. 0)
then 437 do i=1,current_state%number_q_fields
439 q_field_name=trim(
olqbar)//
"_"//trim(q_metadata%name)
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))
450 real_data_1d_double=current_state%global_grid%configuration%vertical%olqbar(:, i))
454 allocate(current_state%global_grid%configuration%vertical%olzqbar(z_size, current_state%number_q_fields))
456 else if (current_state%number_q_fields .gt. 0)
then 457 do i=1,current_state%number_q_fields
459 q_field_name=trim(
olzqbar)//
"_"//trim(q_metadata%name)
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))
470 real_data_1d_double=current_state%global_grid%configuration%vertical%olzqbar(:, i))
482 character(len=*),
intent(in) :: z_key
483 integer,
intent(in) :: z_dim_id, ncid
486 real,
dimension(:),
allocatable :: data
489 allocate(
data(z_size))
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))
499 current_state%global_grid%configuration%vertical%kgd(i) = i
500 current_state%global_grid%configuration%vertical%hgd(i) =
real(data(i))
512 integer,
intent(in) :: ncid, dimension_id, dimension
513 character(len=*),
intent(in) :: variable_key
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))
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))
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))
532 grid%bottom(dimension) = 0
534 grid%size(dimension) = dim_size
535 grid%dimensions = grid%dimensions + 1
536 grid%active(dimension) = .true.
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
559 integer :: variable_id
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 567 if (
present(int_data))
then 572 if (
present(real_data))
then 577 if (
present(real_data_1d))
then 578 if (
present(start) .and.
present(count) .and.
present(map))
then 580 else if (
present(start) .and.
present(count))
then 585 else if (
present(real_data_1d_double))
then 586 if (
present(start) .and.
present(count) .and.
present(map))
then 588 else if (
present(start) .and.
present(count))
then 593 else if (
present(real_data_2d_double))
then 594 if (
present(start) .and.
present(count) .and.
present(map))
then 596 else if (
present(start) .and.
present(count))
then 601 else if (
present(real_data_3d))
then 602 if (
present(start) .and.
present(count) .and.
present(map))
then 604 else if (
present(start) .and.
present(count))
then 609 else if (
present(integer_data_1d))
then 610 if (
present(start) .and.
present(count) .and.
present(map))
then 612 else if (
present(start) .and.
present(count))
then 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
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
integer, parameter, public dual_grid
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.
Contains prognostic field definitions and functions.
character(len= *), parameter zv_key
A prognostic field which is assumed to be 3D.
integer function nf90_inq_dimid(ncid, key, dim_id)
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.
character(len= *), parameter zth_key
integer, parameter, public z_index
Grid index parameters.
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.
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.
subroutine, public log_master_log(level, message)
Will log just from the master process.
integer, parameter nf90_nowrite
Conversion between common inbuilt FORTRAN data types.
integer function nf90_inq_varid(ncid, key, varid)
Converts data types to strings.
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...
Converts data types to logical.
type(q_metadata_type) function, public get_indices_descriptor(i)
Retrieves the indicies descriptor at a specific location.
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...
This manages the Q variables and specifically the mapping between names and the index that they are s...
character(len= *), parameter string_dim_key
String dimension key.
Determines whether a data item can be represented as a logical or not.
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.
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.
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.
subroutine load_misc(current_state, ncid)
Loads in misc data from the checkpoint file.
integer, parameter, public string_length
Default length of strings.
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)
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...
Determines whether a data item can be represented as a real or not.
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.
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.
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.
integer, parameter, public y_index
subroutine, public set_q_index(index, name)
Sets a Q index to be active at a specific index and sets the name.
character(len= *), parameter nqfields
integer, parameter, public x_index
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...
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)