4 use netcdf,
only : nf90_double, nf90_real, nf90_int, nf90_char, nf90_global, nf90_clobber, nf90_netcdf4, nf90_mpiio, &
5 nf90_collective, nf90_def_var, nf90_var_par_access, nf90_def_var_fill, nf90_put_att, nf90_create, nf90_put_var, &
6 nf90_def_dim, nf90_enddef, nf90_close, nf90_inq_dimid, nf90_inq_varid
17 x_dim_key,
y_dim_key,
z_dim_key,
zn_dim_key,
q_dim_key,
q_key,
zq_key,
th_key,
zth_key,
p_key,
u_key,
v_key,
w_key, &
18 zu_key,
zv_key,
zw_key,
x_key,
y_key,
z_key,
zn_key,
nqfields,
ugal,
vgal,
time_key,
timestep,
max_string_length, &
20 q_indices_dim_key,
x_resolution,
y_resolution,
x_top,
y_top,
x_bottom,
y_bottom,
thref,
olubar,
olzubar,
olvbar, &
24 use mpi,
only : mpi_info_null
42 character(len=*),
intent(in) :: filename
44 integer :: ncid, z_dim_id, y_dim_id, x_dim_id, q_dim_id, x_id, y_id, z_id, th_id, p_id, time_id,&
45 u_id, v_id, w_id, q_id, zu_id, zv_id, zw_id, zth_id, zq_id, timestep_id, ugal_id, &
46 vgal_id, number_q_fields_id, string_dim_id, key_value_dim_id, options_id, q_indices_id, &
47 dtm_id, dtm_new_id, absolute_new_dtm_id
48 logical :: q_indices_declared
50 #ifdef SINGLE_MONC_DO_SEQUENTIAL_NETCDF
51 if (current_state%parallel%processes .gt. 1)
then
53 comm = current_state%parallel%monc_communicator, info = mpi_info_null))
59 comm = current_state%parallel%monc_communicator, info = mpi_info_null))
69 if (current_state%number_q_fields .gt. 0)
call define_q_variable(ncid, current_state%parallel%processes .gt. 1, &
70 q_dim_id, z_dim_id, y_dim_id, x_dim_id, q_id, zq_id)
72 x_dim_id, u_id, v_id, w_id, th_id, p_id, zu_id, zv_id, zw_id, zth_id)
74 dtm_id, dtm_new_id, absolute_new_dtm_id)
78 if (current_state%parallel%my_rank==0)
call write_out_grid(ncid, current_state%global_grid)
80 call write_out_all_fields(current_state, ncid, u_id, v_id, w_id, zu_id, zv_id, zw_id, th_id, zth_id, q_id, zq_id, p_id)
81 if (current_state%parallel%my_rank==0)
then
85 ugal_id, vgal_id, number_q_fields_id, dtm_id, dtm_new_id, absolute_new_dtm_id)
90 if (current_state%parallel%my_rank==0)
call write_out_pdf_fields(ncid, current_state%global_grid)
98 integer,
intent(in) :: ncid
100 integer :: date_values(8)
102 call date_and_time(values=date_values)
114 integer,
intent(in) :: ncid, q_indices_id
116 integer :: i, current_index
122 if (specific_q_data%l_used)
then
125 current_index=current_index+1
136 integer,
intent(in) :: ncid, options_id
139 character(len=STRING_LENGTH),
pointer :: sized_raw_character
140 class(*),
pointer :: raw_data, raw_to_string
144 raw_to_string=>raw_data
146 select type (raw_data)
155 type is(
character(len=*))
169 subroutine write_out_all_fields(current_state, ncid, u_id, v_id, w_id, zu_id, zv_id, zw_id, th_id, zth_id, q_id, zq_id, p_id)
171 integer,
intent(in) :: ncid, u_id, v_id, w_id, zu_id, zv_id, zw_id, th_id, zth_id, q_id, zq_id, p_id
174 logical :: multi_process
176 multi_process = current_state%parallel%processes .gt. 1
189 if (current_state%th%active)
then
195 do i=1,current_state%number_q_fields
196 if (current_state%q(i)%active)
then
209 integer,
intent(in) :: ncid, variable_id
212 logical,
intent(in) :: multi_process
213 integer,
optional,
intent(in) :: fourth_dim_loc
215 integer :: start(4), count(4), i, map(4)
217 if (multi_process .or.
present(fourth_dim_loc))
then
222 map(i)=map(i-1)*local_grid%size(i-1)
224 start(i) = local_grid%start(i)
225 count(i) = local_grid%size(i)
227 if (
present(fourth_dim_loc))
then
228 start(4) = fourth_dim_loc
230 map(4)=map(3)*local_grid%size(3)
234 local_grid%local_domain_end_index(
z_index),local_grid%local_domain_start_index(
y_index):&
235 local_grid%local_domain_end_index(
y_index), local_grid%local_domain_start_index(
x_index):&
236 local_grid%local_domain_end_index(
x_index)), start=start, count=count))
239 local_grid%local_domain_end_index(
z_index),local_grid%local_domain_start_index(
y_index):&
240 local_grid%local_domain_end_index(
y_index), local_grid%local_domain_start_index(
x_index):&
241 local_grid%local_domain_end_index(
x_index))))
252 integer,
intent(in) :: ncid
279 integer,
intent(in) :: ncid
284 if (
allocated(grid%configuration%vertical%olubar))
then
288 if (
allocated(grid%configuration%vertical%olzubar))
then
292 if (
allocated(grid%configuration%vertical%olvbar))
then
296 if (
allocated(grid%configuration%vertical%olzvbar))
then
300 if (
allocated(grid%configuration%vertical%olthbar))
then
304 if (
allocated(grid%configuration%vertical%olzthbar))
then
308 if (
allocated(grid%configuration%vertical%olqbar))
then
312 if (
allocated(grid%configuration%vertical%olzqbar))
then
324 integer,
intent(in) :: ncid
326 integer :: z_var_id, zn_var_id, thref_var_id
343 integer,
intent(in) :: ncid
344 integer,
intent(out) :: string_dim_id, key_value_dim_id, options_id
346 integer :: options_dim_id, command_dimensions(3)
352 command_dimensions = (/ string_dim_id, key_value_dim_id, options_dim_id /)
365 integer,
intent(in) :: ncid, string_dim_id, key_value_dim_id
366 integer,
intent(out) :: q_indices_id
368 integer :: q_indices_dim_id, command_dimensions(3), number_active_q
372 if (number_active_q == 0)
then
377 command_dimensions = (/ string_dim_id, key_value_dim_id, q_indices_dim_id /)
390 integer,
intent(in) :: ncid
391 integer,
intent(out) :: q_dim_id
404 integer,
intent(in) :: ncid
405 integer,
intent(out) :: z_dim_id, y_dim_id, x_dim_id
407 integer :: empty_dim_id
411 if (current_state%global_grid%active(
z_index))
then
415 z_dim_id = empty_dim_id
417 if (current_state%global_grid%active(
y_index))
then
420 y_dim_id = empty_dim_id
422 if (current_state%global_grid%active(
x_index))
then
425 x_dim_id = empty_dim_id
440 integer,
intent(in) :: ncid
442 integer :: var_id, z_dim_id
444 if (current_state%global_grid%active(
x_index))
then
452 if (current_state%global_grid%active(
y_index))
then
460 if (current_state%global_grid%active(
z_index))
then
474 integer,
intent(in) :: ncid
476 integer :: var_id, z_dim_id, q_dim_id, qdimids(2)
480 if (
allocated(current_state%global_grid%configuration%vertical%olubar))
then
483 if (
allocated(current_state%global_grid%configuration%vertical%olzubar))
then
486 if (
allocated(current_state%global_grid%configuration%vertical%olvbar))
then
489 if (
allocated(current_state%global_grid%configuration%vertical%olzvbar))
then
492 if (
allocated(current_state%global_grid%configuration%vertical%olthbar))
then
495 if (
allocated(current_state%global_grid%configuration%vertical%olzthbar))
then
498 if (
allocated(current_state%global_grid%configuration%vertical%olqbar) .or. &
499 allocated(current_state%global_grid%configuration%vertical%olzqbar))
then
501 qdimids=(/ z_dim_id, q_dim_id /)
502 if (
allocated(current_state%global_grid%configuration%vertical%olqbar))
then
505 if (
allocated(current_state%global_grid%configuration%vertical%olzqbar))
then
520 subroutine define_q_variable(ncid, multi_process, q_dim_id, z_dim_id, y_dim_id, x_dim_id, q_id, zq_id)
521 logical,
intent(in) :: multi_process
522 integer,
intent(in) :: ncid, z_dim_id, y_dim_id, x_dim_id, q_dim_id
523 integer,
intent(out) :: q_id, zq_id
525 integer,
dimension(:),
allocatable :: dimids
528 dimids = (/ z_dim_id, y_dim_id, x_dim_id, q_dim_id /)
535 if (multi_process)
then
538 call check_status(nf90_var_par_access(ncid, q_id, nf90_collective))
539 call check_status(nf90_var_par_access(ncid, zq_id, nf90_collective))
556 y_dim_id, x_dim_id, u_id, v_id, w_id, th_id, p_id, zu_id, zv_id, zw_id, zth_id)
558 logical,
intent(in) :: multi_process
559 integer,
intent(in) :: ncid, z_dim_id, y_dim_id, x_dim_id
560 integer,
intent(out) :: u_id, v_id, w_id, th_id, p_id, zu_id, zv_id, zw_id, zth_id
574 if (current_state%th%active)
then
578 if (current_state%p%active)
then
587 dtm_id, dtm_new_id, absolute_new_dtm_id)
588 integer,
intent(in) :: ncid
589 integer,
intent(out) :: timestep_id, time_id, ugal_id, vgal_id, number_q_fields_id, dtm_id, dtm_new_id, absolute_new_dtm_id
606 vgal_id, number_q_fields_id, dtm_id, dtm_new_id, absolute_new_dtm_id)
608 integer,
intent(in) :: ncid, timestep_id, time_id, ugal_id, vgal_id, number_q_fields_id, &
609 dtm_id, dtm_new_id, absolute_new_dtm_id
630 integer,
intent(in) :: ncid, dimone
631 integer,
intent(in),
optional :: dimtwo, dimthree
632 integer,
intent(out) :: field_id
633 character(len=*),
intent(in) :: field_name
634 logical,
intent(in) :: multi_process
636 integer,
dimension(:),
allocatable :: dimids
638 if (
present(dimtwo) .and.
present(dimthree))
then
640 dimids = (/ dimone, dimtwo, dimthree /)
641 else if (
present(dimtwo) .or.
present(dimthree))
then
643 dimids = (/ dimone, merge(dimtwo, dimthree,
present(dimtwo)) /)
646 dimids = (/ dimone /)
651 if (multi_process)
then
652 call check_status(nf90_def_var_fill(ncid, field_id, 1, 1))
653 call check_status(nf90_var_par_access(ncid, field_id, nf90_collective))
660 integer,
intent(in) :: ncid
662 integer :: var_id, z_dim_id
666 if (
allocated(current_state%global_grid%configuration%vertical%w_up))
then
669 if (
allocated(current_state%global_grid%configuration%vertical%w_dwn))
then
675 integer,
intent(in) :: ncid
680 if (
allocated(grid%configuration%vertical%w_up))
then
684 if (
allocated(grid%configuration%vertical%w_dwn))
then