MONC
Functions/Subroutines | Variables
stepfields_mod Module Reference

Does the field stepping Stepping is called at the end of processing a column and steps the x-2 column. More...

Functions/Subroutines

type(component_descriptor_type) function, public stepfields_get_descriptor ()
 Provides the descriptor back to the caller and is used in component registration. More...
 
subroutine initialisation_callback (current_state)
 Initialisation callback. More...
 
subroutine finalisation_callback (current_state)
 Finalisation callback. More...
 
subroutine timestep_callback (current_state)
 Called at each timestep and will perform swapping and smoothing as required. More...
 
subroutine step_all_fields (current_state)
 Steps all fields. More...
 
subroutine determine_local_flow_minmax (current_state, local_y, local_x)
 Determines the minimum and maximum values for the local flow field. These are before the stepping, and are all reduced later on in the cfl test. More...
 
subroutine reset_local_minmax_values (current_state)
 Resets the local min and max values for the flow fields. More...
 
subroutine remove_negative_rounding_errors_for_single_field (x_local_index, y_local_index, x_prev, y_prev, field, local_grid)
 Removes the negative rounding errors from a specific single field. This works two columns behind and then catches up on the last column. More...
 
subroutine remove_negative_rounding_errors_in_slice (y_local_index, x_prev, y_prev, field, local_grid)
 Removes the negative rounding errors from a slice of a single field. This works two columns behind and then catches up on the last column. More...
 
subroutine step_single_field (x_local_index, y_local_index, x_prev, y_prev, field, zfield, sfield, local_grid, flow_field, direction, dtm, gal, c1, c2, do_timesmoothing, sav)
 Steps a single specific field. This will step on the yth column of the x-2 slice and x-1 and x if this is the last slice. More...
 
subroutine perform_timesmooth_for_field (field, zfield, local_grid, x_index, y_index, c1, c2)
 Performs initial timesmoothing for a theta or Q field using Robert filter. This is finished off in swapsmooth. More...
 
subroutine step_column_in_slice (y_local_index, x_prev, y_prev, field, zfield, sfield, local_grid, flow_field, direction, dtm, gal, c1, c2, do_timesmoothing, sav)
 Will step a column in a specific slice. If y_prev is large enough then will step the y-1 column and if this is the last column of the slice then will also step the current column. More...
 
subroutine step_field (x_local_index, y_local_index, field, zfield, sfield, local_grid, flow_field, direction, dtm, gal, sav)
 Will do the actual field stepping. More...
 
subroutine compute_component_tendencies (current_state, cxn, cyn, txn, tyn)
 Computation of component tendencies. More...
 
subroutine field_information_retrieval_callback (current_state, name, field_information)
 Field information retrieval callback, this returns information for a specific component's published field. More...
 
subroutine field_value_retrieval_callback (current_state, name, field_value)
 Field value retrieval callback, this returns the value of a specific published field. More...
 
subroutine set_published_field_value (field_value, real_1d_field, real_2d_field, real_3d_field)
 Sets the published field value from the temporary diagnostic values held by this component. More...
 

Variables

logical determine_flow_minmax =.false.
 
logical cfl_is_enabled
 
real(kind=default_precision), dimension(:), allocatable resetq_min
 
logical l_nonconservative_positive_q =.true.
 
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_th
 
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_qv
 
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_ql
 
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_qi
 
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_qr
 
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_qs
 
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_qg
 
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_tabs
 
logical l_tend_3d_th
 
logical l_tend_3d_qv
 
logical l_tend_3d_ql
 
logical l_tend_3d_qi
 
logical l_tend_3d_qr
 
logical l_tend_3d_qs
 
logical l_tend_3d_qg
 
logical l_tend_3d_tabs
 
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_th
 
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_qv
 
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_ql
 
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_qi
 
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_qr
 
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_qs
 
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_qg
 
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_tabs
 
logical l_tend_pr_tot_th
 
logical l_tend_pr_tot_qv
 
logical l_tend_pr_tot_ql
 
logical l_tend_pr_tot_qi
 
logical l_tend_pr_tot_qr
 
logical l_tend_pr_tot_qs
 
logical l_tend_pr_tot_qg
 
logical l_tend_pr_tot_tabs
 
integer iqv =0
 
integer iql =0
 
integer iqr =0
 
integer iqi =0
 
integer iqs =0
 
integer iqg =0
 
integer diagnostic_generation_frequency
 

Detailed Description

Does the field stepping Stepping is called at the end of processing a column and steps the x-2 column.

Function/Subroutine Documentation

◆ compute_component_tendencies()

subroutine stepfields_mod::compute_component_tendencies ( type(model_state_type), intent(inout), target  current_state,
integer, intent(in)  cxn,
integer, intent(in)  cyn,
integer, intent(in)  txn,
integer, intent(in)  tyn 
)
private

Computation of component tendencies.

Parameters
current_stateCurrent model state
cxnThe current slice, x, index
cynThe current column, y, index.
txntarget_x_index
tyntarget_y_index

Definition at line 628 of file stepfields.F90.

629  type(model_state_type), target, intent(inout) :: current_state
630  integer, intent(in) :: cxn, cyn, txn, tyn
631 
632  ! Calculate change in tendency due to component
633  if (l_tend_3d_th) then
634  tend_3d_th(:,tyn,txn)=current_state%sth%data(:,cyn,cxn)
635  endif
636  if (l_tend_3d_qv) then
637  tend_3d_qv(:,tyn,txn)=current_state%sq(iqv)%data(:,cyn,cxn)
638  endif
639  if (l_tend_3d_ql) then
640  tend_3d_ql(:,tyn,txn)=current_state%sq(iql)%data(:,cyn,cxn)
641  endif
642  if (l_tend_3d_qi) then
643  tend_3d_qi(:,tyn,txn)=current_state%sq(iqi)%data(:,cyn,cxn)
644  endif
645  if (l_tend_3d_qr) then
646  tend_3d_qr(:,tyn,txn)=current_state%sq(iqr)%data(:,cyn,cxn)
647  endif
648  if (l_tend_3d_qs) then
649  tend_3d_qs(:,tyn,txn)=current_state%sq(iqs)%data(:,cyn,cxn)
650  endif
651  if (l_tend_3d_qg) then
652  tend_3d_qg(:,tyn,txn)=current_state%sq(iqg)%data(:,cyn,cxn)
653  endif
654  if (l_tend_3d_tabs) then
655  tend_3d_tabs(:,tyn,txn)= &
656  current_state%sth%data(:,cyn,cxn) * current_state%global_grid%configuration%vertical%rprefrcp(:)
657  endif
658 
659  ! Add local tendency fields to the profile total
660  if (l_tend_pr_tot_th) then
661  tend_pr_tot_th(:)=tend_pr_tot_th(:) + tend_3d_th(:,tyn,txn)
662  endif
663  if (l_tend_pr_tot_qv) then
664  tend_pr_tot_qv(:)=tend_pr_tot_qv(:) + tend_3d_qv(:,tyn,txn)
665  endif
666  if (l_tend_pr_tot_ql) then
667  tend_pr_tot_ql(:)=tend_pr_tot_ql(:) + tend_3d_ql(:,tyn,txn)
668  endif
669  if (l_tend_pr_tot_qi) then
670  tend_pr_tot_qi(:)=tend_pr_tot_qi(:) + tend_3d_qi(:,tyn,txn)
671  endif
672  if (l_tend_pr_tot_qr) then
673  tend_pr_tot_qr(:)=tend_pr_tot_qr(:) + tend_3d_qr(:,tyn,txn)
674  endif
675  if (l_tend_pr_tot_qs) then
676  tend_pr_tot_qs(:)=tend_pr_tot_qs(:) + tend_3d_qs(:,tyn,txn)
677  endif
678  if (l_tend_pr_tot_qg) then
679  tend_pr_tot_qg(:)=tend_pr_tot_qg(:) + tend_3d_qg(:,tyn,txn)
680  endif
681  if (l_tend_pr_tot_tabs) then
682  tend_pr_tot_tabs(:)=tend_pr_tot_tabs(:) + tend_3d_tabs(:,tyn,txn)
683  endif
684 
Here is the caller graph for this function:

◆ determine_local_flow_minmax()

subroutine stepfields_mod::determine_local_flow_minmax ( type(model_state_type), intent(inout), target  current_state,
integer, intent(in)  local_y,
integer, intent(in)  local_x 
)
private

Determines the minimum and maximum values for the local flow field. These are before the stepping, and are all reduced later on in the cfl test.

Parameters
current_stateThe current model state
local_yThe local y index
local_xThe local x index

Definition at line 365 of file stepfields.F90.

366  type(model_state_type), target, intent(inout) :: current_state
367  integer, intent(in) :: local_x, local_y
368 
369  integer :: k
370 
371  do k=2, current_state%local_grid%local_domain_end_index(z_index)
372 #ifdef U_ACTIVE
373  current_state%local_zumax = max(current_state%local_zumax, current_state%zu%data(k,local_y,local_x))
374  current_state%local_zumin = min(current_state%local_zumin, current_state%zu%data(k,local_y,local_x))
375 #endif
376 #ifdef V_ACTIVE
377  current_state%local_zvmax = max(current_state%local_zvmax, current_state%zv%data(k,local_y,local_x))
378  current_state%local_zvmin = min(current_state%local_zvmin, current_state%zv%data(k,local_y,local_x))
379 #endif
380 #ifdef W_ACTIVE
381  if (k .lt. current_state%local_grid%local_domain_end_index(z_index)) then
382  current_state%abswmax(k) = max(current_state%abswmax(k), abs(current_state%zw%data(k,local_y,local_x)))
383  end if
384 #endif
385  end do
Here is the caller graph for this function:

◆ field_information_retrieval_callback()

subroutine stepfields_mod::field_information_retrieval_callback ( type(model_state_type), intent(inout), target  current_state,
character(len=*), intent(in)  name,
type(component_field_information_type), intent(out)  field_information 
)
private

Field information retrieval callback, this returns information for a specific component's published field.

Parameters
current_stateCurrent model state
nameThe name of the field to retrieve information for
field_informationPopulated with information about the field
strcompStarting index within 1st argument string that matches substring (2nd argument); 0 if not a match.

Definition at line 693 of file stepfields.F90.

694  type(model_state_type), target, intent(inout) :: current_state
695  character(len=*), intent(in) :: name
696  type(component_field_information_type), intent(out) :: field_information
697  integer :: strcomp
698 
699  ! Field information for 3d
700  strcomp=index(name, "_total_3d_local")
701  if (strcomp .ne. 0) then
702  field_information%field_type=component_array_field_type
703  field_information%number_dimensions=3
704  field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
705  field_information%dimension_sizes(2)=current_state%local_grid%size(y_index)
706  field_information%dimension_sizes(3)=current_state%local_grid%size(x_index)
707  field_information%data_type=component_double_data_type
708 
709  if (name .eq. "tend_th_total_3d_local") then
710  field_information%enabled=l_tend_3d_th
711  else if (name .eq. "tend_qv_total_3d_local") then
712  field_information%enabled=l_tend_3d_qv
713  else if (name .eq. "tend_ql_total_3d_local") then
714  field_information%enabled=l_tend_3d_ql
715  else if (name .eq. "tend_qi_total_3d_local") then
716  field_information%enabled=l_tend_3d_qi
717  else if (name .eq. "tend_qr_total_3d_local") then
718  field_information%enabled=l_tend_3d_qr
719  else if (name .eq. "tend_qs_total_3d_local") then
720  field_information%enabled=l_tend_3d_qs
721  else if (name .eq. "tend_qg_total_3d_local") then
722  field_information%enabled=l_tend_3d_qg
723  else if (name .eq. "tend_tabs_total_3d_local") then
724  field_information%enabled=l_tend_3d_tabs
725  else
726  field_information%enabled=.true.
727  end if
728 
729  end if !end 3d check
730 
731  ! Field information for profiles
732  strcomp=index(name, "_total_profile_total_local")
733  if (strcomp .ne. 0) then
734  field_information%field_type=component_array_field_type
735  field_information%number_dimensions=1
736  field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
737  field_information%data_type=component_double_data_type
738 
739  if (name .eq. "tend_th_total_profile_total_local") then
740  field_information%enabled=l_tend_pr_tot_th
741  else if (name .eq. "tend_qv_total_profile_total_local") then
742  field_information%enabled=l_tend_pr_tot_qv
743  else if (name .eq. "tend_ql_total_profile_total_local") then
744  field_information%enabled=l_tend_pr_tot_ql
745  else if (name .eq. "tend_qi_total_profile_total_local") then
746  field_information%enabled=l_tend_pr_tot_qi
747  else if (name .eq. "tend_qr_total_profile_total_local") then
748  field_information%enabled=l_tend_pr_tot_qr
749  else if (name .eq. "tend_qs_total_profile_total_local") then
750  field_information%enabled=l_tend_pr_tot_qs
751  else if (name .eq. "tend_qg_total_profile_total_local") then
752  field_information%enabled=l_tend_pr_tot_qg
753  else if (name .eq. "tend_tabs_total_profile_total_local") then
754  field_information%enabled=l_tend_pr_tot_tabs
755  else
756  field_information%enabled=.true.
757  end if
758 
759  end if !end profile check
760 
Here is the caller graph for this function:

◆ field_value_retrieval_callback()

subroutine stepfields_mod::field_value_retrieval_callback ( type(model_state_type), intent(inout), target  current_state,
character(len=*), intent(in)  name,
type(component_field_value_type), intent(out)  field_value 
)
private

Field value retrieval callback, this returns the value of a specific published field.

Parameters
current_stateCurrent model state
nameThe name of the field to retrieve the value for
field_valuePopulated with the value of the field

Definition at line 768 of file stepfields.F90.

769  type(model_state_type), target, intent(inout) :: current_state
770  character(len=*), intent(in) :: name
771  type(component_field_value_type), intent(out) :: field_value
772 
773  ! 3d Tendency Fields
774  if (name .eq. "tend_th_total_3d_local" .and. allocated(tend_3d_th)) then
775  call set_published_field_value(field_value, real_3d_field=tend_3d_th)
776  else if (name .eq. "tend_qv_total_3d_local" .and. allocated(tend_3d_qv)) then
777  call set_published_field_value(field_value, real_3d_field=tend_3d_qv)
778  else if (name .eq. "tend_ql_total_3d_local" .and. allocated(tend_3d_ql)) then
779  call set_published_field_value(field_value, real_3d_field=tend_3d_ql)
780  else if (name .eq. "tend_qi_total_3d_local" .and. allocated(tend_3d_qi)) then
781  call set_published_field_value(field_value, real_3d_field=tend_3d_qi)
782  else if (name .eq. "tend_qr_total_3d_local" .and. allocated(tend_3d_qr)) then
783  call set_published_field_value(field_value, real_3d_field=tend_3d_qr)
784  else if (name .eq. "tend_qs_total_3d_local" .and. allocated(tend_3d_qs)) then
785  call set_published_field_value(field_value, real_3d_field=tend_3d_qs)
786  else if (name .eq. "tend_qg_total_3d_local" .and. allocated(tend_3d_qg)) then
787  call set_published_field_value(field_value, real_3d_field=tend_3d_qg)
788  else if (name .eq. "tend_tabs_total_3d_local" .and. allocated(tend_3d_tabs)) then
789  call set_published_field_value(field_value, real_3d_field=tend_3d_tabs)
790 
791  ! Profile Tendency Fields
792  else if (name .eq. "tend_th_total_profile_total_local" .and. allocated(tend_pr_tot_th)) then
793  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_th)
794  else if (name .eq. "tend_qv_total_profile_total_local" .and. allocated(tend_pr_tot_qv)) then
795  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_qv)
796  else if (name .eq. "tend_ql_total_profile_total_local" .and. allocated(tend_pr_tot_ql)) then
797  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_ql)
798  else if (name .eq. "tend_qi_total_profile_total_local" .and. allocated(tend_pr_tot_qi)) then
799  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_qi)
800  else if (name .eq. "tend_qr_total_profile_total_local" .and. allocated(tend_pr_tot_qr)) then
801  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_qr)
802  else if (name .eq. "tend_qs_total_profile_total_local" .and. allocated(tend_pr_tot_qs)) then
803  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_qs)
804  else if (name .eq. "tend_qg_total_profile_total_local" .and. allocated(tend_pr_tot_qg)) then
805  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_qg)
806  else if (name .eq. "tend_tabs_total_profile_total_local" .and. allocated(tend_pr_tot_tabs)) then
807  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_tabs)
808  end if
809 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ finalisation_callback()

subroutine stepfields_mod::finalisation_callback ( type(model_state_type), intent(inout), target  current_state)
private

Finalisation callback.

Parameters
current_stateThe current model state

Definition at line 209 of file stepfields.F90.

210  type(model_state_type), target, intent(inout) :: current_state
211 
212  deallocate(resetq_min)
213 
214  if (allocated(tend_3d_th)) deallocate(tend_3d_th)
215  if (allocated(tend_3d_qv)) deallocate(tend_3d_qv)
216  if (allocated(tend_3d_ql)) deallocate(tend_3d_ql)
217  if (allocated(tend_3d_qi)) deallocate(tend_3d_qi)
218  if (allocated(tend_3d_qr)) deallocate(tend_3d_qr)
219  if (allocated(tend_3d_qs)) deallocate(tend_3d_qs)
220  if (allocated(tend_3d_qg)) deallocate(tend_3d_qg)
221  if (allocated(tend_3d_tabs)) deallocate(tend_3d_tabs)
222 
223  if (allocated(tend_pr_tot_th)) deallocate(tend_pr_tot_th)
224  if (allocated(tend_pr_tot_qv)) deallocate(tend_pr_tot_qv)
225  if (allocated(tend_pr_tot_ql)) deallocate(tend_pr_tot_ql)
226  if (allocated(tend_pr_tot_qi)) deallocate(tend_pr_tot_qi)
227  if (allocated(tend_pr_tot_qr)) deallocate(tend_pr_tot_qr)
228  if (allocated(tend_pr_tot_qs)) deallocate(tend_pr_tot_qs)
229  if (allocated(tend_pr_tot_qg)) deallocate(tend_pr_tot_qg)
230  if (allocated(tend_pr_tot_tabs)) deallocate(tend_pr_tot_tabs)
231 
Here is the caller graph for this function:

◆ initialisation_callback()

subroutine stepfields_mod::initialisation_callback ( type(model_state_type), intent(inout), target  current_state)
private

Initialisation callback.

Parameters
current_stateThe current model state

Definition at line 89 of file stepfields.F90.

90  type(model_state_type), target, intent(inout) :: current_state
91  logical :: l_qdiag
92 
93 #ifdef U_ACTIVE
94  allocate(current_state%global_grid%configuration%vertical%savolubar(current_state%local_grid%size(z_index)))
95 #endif
96 #ifdef V_ACTIVE
97  allocate(current_state%global_grid%configuration%vertical%savolvbar(current_state%local_grid%size(z_index)))
98 #endif
99 
100  allocate(resetq_min(current_state%number_q_fields))
101  cfl_is_enabled=is_component_enabled(current_state%options_database, "cfltest")
102  if (cfl_is_enabled) call reset_local_minmax_values(current_state)
103 
104  ! Set tendency diagnostic logicals based on availability
105  ! Need to use 3d tendencies to compute the profiles, so they will be allocated
106  ! in the case where profiles are available
107  l_qdiag = (.not. current_state%passive_q .and. current_state%number_q_fields .gt. 0)
108 
109  l_tend_pr_tot_th = current_state%th%active
110  l_tend_pr_tot_qv = l_qdiag .and. current_state%number_q_fields .ge. 1
111  l_tend_pr_tot_ql = l_qdiag .and. current_state%number_q_fields .ge. 2
112  l_tend_pr_tot_qi = l_qdiag .and. current_state%number_q_fields .ge. 11
113  l_tend_pr_tot_qr = l_qdiag .and. current_state%number_q_fields .ge. 11
114  l_tend_pr_tot_qs = l_qdiag .and. current_state%number_q_fields .ge. 11
115  l_tend_pr_tot_qg = l_qdiag .and. current_state%number_q_fields .ge. 11
116  l_tend_pr_tot_tabs = l_tend_pr_tot_th
117 
118  l_tend_3d_th = current_state%th%active .or. l_tend_pr_tot_th
119  l_tend_3d_qv = (l_qdiag .and. current_state%number_q_fields .ge. 1) .or. l_tend_pr_tot_qv
120  l_tend_3d_ql = (l_qdiag .and. current_state%number_q_fields .ge. 2) .or. l_tend_pr_tot_ql
121  l_tend_3d_qi = (l_qdiag .and. current_state%number_q_fields .ge. 11) .or. l_tend_pr_tot_qi
122  l_tend_3d_qr = (l_qdiag .and. current_state%number_q_fields .ge. 11) .or. l_tend_pr_tot_qr
123  l_tend_3d_qs = (l_qdiag .and. current_state%number_q_fields .ge. 11) .or. l_tend_pr_tot_qs
124  l_tend_3d_qg = (l_qdiag .and. current_state%number_q_fields .ge. 11) .or. l_tend_pr_tot_qg
125  l_tend_3d_tabs = l_tend_3d_th
126 
127  ! Allocate 3d tendency fields upon availability
128  if (l_tend_3d_th) then
129  allocate( tend_3d_th(current_state%local_grid%size(z_index), &
130  current_state%local_grid%size(y_index), &
131  current_state%local_grid%size(x_index) ) )
132  endif
133  if (l_tend_3d_qv) then
134  iqv=get_q_index(standard_q_names%VAPOUR, 'stepfields')
135  allocate( tend_3d_qv(current_state%local_grid%size(z_index), &
136  current_state%local_grid%size(y_index), &
137  current_state%local_grid%size(x_index) ) )
138  endif
139  if (l_tend_3d_ql) then
140  iql=get_q_index(standard_q_names%CLOUD_LIQUID_MASS, 'stepfields')
141  allocate( tend_3d_ql(current_state%local_grid%size(z_index), &
142  current_state%local_grid%size(y_index), &
143  current_state%local_grid%size(x_index) ) )
144  endif
145  if (l_tend_3d_qi) then
146  iqi=get_q_index(standard_q_names%ICE_MASS, 'stepfields')
147  allocate( tend_3d_qi(current_state%local_grid%size(z_index), &
148  current_state%local_grid%size(y_index), &
149  current_state%local_grid%size(x_index) ) )
150  endif
151  if (l_tend_3d_qr) then
152  iqr=get_q_index(standard_q_names%RAIN_MASS, 'stepfields')
153  allocate( tend_3d_qr(current_state%local_grid%size(z_index), &
154  current_state%local_grid%size(y_index), &
155  current_state%local_grid%size(x_index) ) )
156  endif
157  if (l_tend_3d_qs) then
158  iqs=get_q_index(standard_q_names%SNOW_MASS, 'stepfields')
159  allocate( tend_3d_qs(current_state%local_grid%size(z_index), &
160  current_state%local_grid%size(y_index), &
161  current_state%local_grid%size(x_index) ) )
162  endif
163  if (l_tend_3d_qg) then
164  iqg=get_q_index(standard_q_names%GRAUPEL_MASS, 'stepfields')
165  allocate( tend_3d_qg(current_state%local_grid%size(z_index), &
166  current_state%local_grid%size(y_index), &
167  current_state%local_grid%size(x_index) ) )
168  endif
169  if (l_tend_3d_tabs) then
170  allocate( tend_3d_tabs(current_state%local_grid%size(z_index), &
171  current_state%local_grid%size(y_index), &
172  current_state%local_grid%size(x_index) ) )
173  endif
174 
175  ! Allocate profile tendency fields upon availability
176  if (l_tend_pr_tot_th) then
177  allocate( tend_pr_tot_th(current_state%local_grid%size(z_index)) )
178  endif
179  if (l_tend_pr_tot_qv) then
180  allocate( tend_pr_tot_qv(current_state%local_grid%size(z_index)) )
181  endif
182  if (l_tend_pr_tot_ql) then
183  allocate( tend_pr_tot_ql(current_state%local_grid%size(z_index)) )
184  endif
185  if (l_tend_pr_tot_qi) then
186  allocate( tend_pr_tot_qi(current_state%local_grid%size(z_index)) )
187  endif
188  if (l_tend_pr_tot_qr) then
189  allocate( tend_pr_tot_qr(current_state%local_grid%size(z_index)) )
190  endif
191  if (l_tend_pr_tot_qs) then
192  allocate( tend_pr_tot_qs(current_state%local_grid%size(z_index)) )
193  endif
194  if (l_tend_pr_tot_qg) then
195  allocate( tend_pr_tot_qg(current_state%local_grid%size(z_index)) )
196  endif
197  if (l_tend_pr_tot_tabs) then
198  allocate( tend_pr_tot_tabs(current_state%local_grid%size(z_index)) )
199  endif
200 
201  ! Save the sampling_frequency to force diagnostic calculation on select time steps
202  diagnostic_generation_frequency=options_get_integer(current_state%options_database, "sampling_frequency")
203 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ perform_timesmooth_for_field()

subroutine stepfields_mod::perform_timesmooth_for_field ( type(prognostic_field_type), intent(inout)  field,
type(prognostic_field_type), intent(inout)  zfield,
type(local_grid_type), intent(inout)  local_grid,
integer, intent(in)  x_index,
integer, intent(in)  y_index,
real(kind=default_precision), intent(in)  c1,
real(kind=default_precision), intent(in)  c2 
)
private

Performs initial timesmoothing for a theta or Q field using Robert filter. This is finished off in swapsmooth.

Parameters
fieldThe field to smooth
zfieldThe zfield to use in smoothing
local_gridDescription of the local grid
x_indexThe X index to work on
y_indexThe Y index to work on
c1Constant to use in smoothing
c2Constant to use in smoothing

Definition at line 515 of file stepfields.F90.

516  type(prognostic_field_type), intent(inout) :: field, zfield
517  type(local_grid_type), intent(inout) :: local_grid
518  integer, intent(in) :: x_index, y_index
519  real(kind=default_precision), intent(in) :: c1, c2
520 
521  integer :: k
522 
523  do k=1,local_grid%size(z_index)
524  field%data(k, y_index, x_index)=c1*field%data(k, y_index, x_index)+c2*zfield%data(k, y_index, x_index)
525  end do
Here is the caller graph for this function:

◆ remove_negative_rounding_errors_for_single_field()

subroutine stepfields_mod::remove_negative_rounding_errors_for_single_field ( integer, intent(in)  x_local_index,
integer, intent(in)  y_local_index,
integer, intent(in)  x_prev,
integer, intent(in)  y_prev,
type(prognostic_field_type), intent(inout)  field,
type(local_grid_type), intent(inout)  local_grid 
)
private

Removes the negative rounding errors from a specific single field. This works two columns behind and then catches up on the last column.

Parameters
x_local_indexThe current local x index
y_local_indexThe current local y index
x_prevThe previous x index to step
y_prevThe previous y index to step
fieldThe prognostic field
local_gridDescription of the local grid

Definition at line 409 of file stepfields.F90.

410  integer, intent(in) :: x_local_index, y_local_index, x_prev, y_prev
411  type(local_grid_type), intent(inout) :: local_grid
412  type(prognostic_field_type), intent(inout) :: field
413 
414  if (x_prev .ge. local_grid%local_domain_start_index(x_index)) then
415  call remove_negative_rounding_errors_in_slice(y_local_index, x_prev, y_prev, field, local_grid)
416  end if
417  if (x_local_index == local_grid%local_domain_end_index(x_index)) then
418  if (x_local_index .gt. 1) then
419  call remove_negative_rounding_errors_in_slice(y_local_index, x_local_index-1, y_prev, field, local_grid)
420  end if
421  call remove_negative_rounding_errors_in_slice(y_local_index, x_local_index, y_prev, field, local_grid)
422  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ remove_negative_rounding_errors_in_slice()

subroutine stepfields_mod::remove_negative_rounding_errors_in_slice ( integer, intent(in)  y_local_index,
integer, intent(in)  x_prev,
integer, intent(in)  y_prev,
type(prognostic_field_type), intent(inout)  field,
type(local_grid_type), intent(inout)  local_grid 
)
private

Removes the negative rounding errors from a slice of a single field. This works two columns behind and then catches up on the last column.

Parameters
y_local_indexThe current local y index
x_prevThe previous x index to step
y_prevThe previous y index to step
fieldThe prognostic field
local_gridDescription of the local grid

Definition at line 432 of file stepfields.F90.

433  integer, intent(in) :: y_local_index, x_prev, y_prev
434  type(local_grid_type), intent(inout) :: local_grid
435  type(prognostic_field_type), intent(inout) :: field
436 
437  if (y_prev .ge. local_grid%local_domain_start_index(y_index)) then
438  where (field%data(:, y_prev, x_prev) < 0.0_default_precision)
439  field%data(:, y_prev, x_prev)=0.0_default_precision
440  end where
441  end if
442  if (y_local_index == local_grid%local_domain_end_index(y_index)) then
443  where (field%data(:, y_local_index, x_prev) < 0.0_default_precision)
444  field%data(:, y_local_index, x_prev)=0.0_default_precision
445  end where
446  end if
Here is the caller graph for this function:

◆ reset_local_minmax_values()

subroutine stepfields_mod::reset_local_minmax_values ( type(model_state_type), intent(inout), target  current_state)
private

Resets the local min and max values for the flow fields.

Parameters
current_stateThe current model state

Definition at line 390 of file stepfields.F90.

391  type(model_state_type), intent(inout), target :: current_state
392 
393  ! Reset the local values for the next timestep
394  current_state%local_zumin=rlargep
395  current_state%local_zumax=-rlargep
396  current_state%local_zvmin=rlargep
397  current_state%local_zvmax=-rlargep
398  current_state%abswmax=-rlargep
Here is the caller graph for this function:

◆ set_published_field_value()

subroutine stepfields_mod::set_published_field_value ( type(component_field_value_type), intent(inout)  field_value,
real(kind=default_precision), dimension(:), optional  real_1d_field,
real(kind=default_precision), dimension(:,:), optional  real_2d_field,
real(kind=default_precision), dimension(:,:,:), optional  real_3d_field 
)
private

Sets the published field value from the temporary diagnostic values held by this component.

Parameters
field_valuePopulated with the value of the field
real_1d_fieldOptional one dimensional real of values to publish
real_2d_fieldOptional two dimensional real of values to publish

Definition at line 817 of file stepfields.F90.

818  type(component_field_value_type), intent(inout) :: field_value
819  real(kind=default_precision), dimension(:), optional :: real_1d_field
820  real(kind=default_precision), dimension(:,:), optional :: real_2d_field
821  real(kind=default_precision), dimension(:,:,:), optional :: real_3d_field
822 
823  if (present(real_1d_field)) then
824  allocate(field_value%real_1d_array(size(real_1d_field)), source=real_1d_field)
825  else if (present(real_2d_field)) then
826  allocate(field_value%real_2d_array(size(real_2d_field, 1), size(real_2d_field, 2)), source=real_2d_field)
827  else if (present(real_3d_field)) then
828  allocate(field_value%real_3d_array(size(real_3d_field, 1), size(real_3d_field, 2), size(real_3d_field, 3)), &
829  source=real_3d_field)
830  end if
Here is the caller graph for this function:

◆ step_all_fields()

subroutine stepfields_mod::step_all_fields ( type(model_state_type), intent(inout), target  current_state)
private

Steps all fields.

Parameters
current_stateThe current model state_mod

Definition at line 315 of file stepfields.F90.

316  type(model_state_type), target, intent(inout) :: current_state
317 
318  integer :: x_prev, y_prev, i, k
319  real(kind=default_precision) :: c1, c2
320 
321  x_prev = current_state%column_local_x-2
322  y_prev = current_state%column_local_y-1
323 
324  c1 = 1.0_default_precision - 2.0_default_precision*current_state%tsmth
325  c2 = current_state%tsmth
326 
327 #ifdef U_ACTIVE
328  current_state%global_grid%configuration%vertical%savolubar=current_state%global_grid%configuration%vertical%olzubar
329  call step_single_field(current_state%column_local_x, current_state%column_local_y, &
330  x_prev, y_prev, current_state%u, current_state%zu, current_state%su, current_state%local_grid, .true., &
331  current_state%field_stepping, current_state%dtm, current_state%ugal, c1, c2, .false., current_state%savu)
332 #endif
333 #ifdef V_ACTIVE
334  current_state%global_grid%configuration%vertical%savolvbar=current_state%global_grid%configuration%vertical%olzvbar
335  call step_single_field(current_state%column_local_x, current_state%column_local_y, &
336  x_prev, y_prev, current_state%v, current_state%zv, current_state%sv, current_state%local_grid, .true., &
337  current_state%field_stepping, current_state%dtm, current_state%vgal, c1, c2, .false., current_state%savv)
338 #endif
339 #ifdef W_ACTIVE
340  call step_single_field(current_state%column_local_x, current_state%column_local_y, &
341  x_prev, y_prev, current_state%w, current_state%zw, current_state%sw, current_state%local_grid, .false., &
342  current_state%field_stepping, current_state%dtm, real(0., kind=default_precision), c1, c2, .false., current_state%savw)
343 #endif
344  if (current_state%th%active) then
345  call step_single_field(current_state%column_local_x, current_state%column_local_y, &
346  x_prev, y_prev, current_state%th, current_state%zth, current_state%sth, current_state%local_grid, .false., &
347  current_state%field_stepping, current_state%dtm, real(0., kind=default_precision), c1, c2, &
348  current_state%field_stepping == centred_stepping)
349  endif
350  do i=1,current_state%number_q_fields
351  if (current_state%q(i)%active) then
352  call step_single_field(current_state%column_local_x, current_state%column_local_y, x_prev, y_prev, &
353  current_state%q(i), current_state%zq(i), current_state%sq(i), current_state%local_grid, .false., &
354  current_state%field_stepping, current_state%dtm, real(0., kind=default_precision), c1, c2, &
355  current_state%field_stepping == centred_stepping)
356  end if
357  end do
Here is the call graph for this function:
Here is the caller graph for this function:

◆ step_column_in_slice()

subroutine stepfields_mod::step_column_in_slice ( integer, intent(in)  y_local_index,
integer, intent(in)  x_prev,
integer, intent(in)  y_prev,
type(prognostic_field_type), intent(inout)  field,
type(prognostic_field_type), intent(inout)  zfield,
type(prognostic_field_type), intent(inout)  sfield,
type(local_grid_type), intent(inout)  local_grid,
logical, intent(in)  flow_field,
integer, intent(in)  direction,
real(kind=default_precision), intent(in)  dtm,
real(kind=default_precision), intent(in)  gal,
real(kind=default_precision), intent(in)  c1,
real(kind=default_precision), intent(in)  c2,
logical, intent(in)  do_timesmoothing,
type(prognostic_field_type), intent(inout), optional  sav 
)
private

Will step a column in a specific slice. If y_prev is large enough then will step the y-1 column and if this is the last column of the slice then will also step the current column.

Parameters
x_local_indexThe current local x index
y_local_indexThe current local y index
x_prevThe previous x index to step
y_prevThe previous y index to step
fieldThe prognostic field
zfieldZ prognostic field
sfieldSource terms for the prognostic field
local_gridDescription of the local grid
flow_fieldWhether or not this is a flow field
directionThe stepping direction (centred or forward)
dtmThe delta time per timestep
galGalilean transformation
c1Constant to use in smoothing
c2Constant to use in smoothing
do_timesmoothingWhether timesmoothing using Robert filter should be done on the field
savOptional sav field

Definition at line 546 of file stepfields.F90.

548  integer, intent(in) :: y_local_index, x_prev, y_prev, direction
549  real(kind=default_precision), intent(in) :: dtm, gal
550  logical, intent(in) :: flow_field, do_timesmoothing
551  type(local_grid_type), intent(inout) :: local_grid
552  type(prognostic_field_type), intent(inout) :: field, zfield, sfield
553  real(kind=default_precision), intent(in) :: c1, c2
554  type(prognostic_field_type), optional, intent(inout) :: sav
555 
556  if (y_prev .ge. local_grid%local_domain_start_index(y_index)) then
557  if (do_timesmoothing) then
558  call perform_timesmooth_for_field(field, zfield, local_grid, x_prev, y_prev, c1, c2)
559  end if
560  if (present(sav)) then
561  call step_field(x_prev, y_prev, field, zfield, sfield, local_grid, flow_field, direction, dtm, gal, sav)
562  else
563  call step_field(x_prev, y_prev, field, zfield, sfield, local_grid, flow_field, direction, dtm, gal)
564  end if
565  end if
566 
567  if (y_local_index == local_grid%local_domain_end_index(y_index)) then
568  if (do_timesmoothing) then
569  call perform_timesmooth_for_field(field, zfield, local_grid, x_prev, y_local_index, c1, c2)
570  end if
571  if (present(sav)) then
572  call step_field(x_prev, y_local_index, field, zfield, sfield, local_grid, flow_field, direction, dtm, gal, sav)
573  else
574  call step_field(x_prev, y_local_index, field, zfield, sfield, local_grid, flow_field, direction, dtm, gal)
575  end if
576  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ step_field()

subroutine stepfields_mod::step_field ( integer, intent(in)  x_local_index,
integer, intent(in)  y_local_index,
type(prognostic_field_type), intent(inout)  field,
type(prognostic_field_type), intent(inout)  zfield,
type(prognostic_field_type), intent(inout)  sfield,
type(local_grid_type), intent(inout)  local_grid,
logical, intent(in)  flow_field,
integer, intent(in)  direction,
real(kind=default_precision), intent(in)  dtm,
real(kind=default_precision), intent(in)  gal,
type(prognostic_field_type), intent(inout), optional  sav 
)
private

Will do the actual field stepping.

Parameters
flow_fieldWhether or not we are stepping a flow field
direction1=forward, 0=centred
x_indexThe local X slice index
y_indexThe local Y column index
kkpPoints in the vertical column
dtmThe model timestep
fieldThe prognostic field
zfieldThe prognostic z field (which goes to timestep t+1)
xfieldThe tendency of the field
galThe galilean transformation

Definition at line 590 of file stepfields.F90.

591  integer, intent(in) :: x_local_index, y_local_index, direction
592  real(kind=default_precision), intent(in) :: dtm, gal
593  logical, intent(in) :: flow_field
594  type(local_grid_type), intent(inout) :: local_grid
595  type(prognostic_field_type), intent(inout) :: field, zfield, sfield
596  type(prognostic_field_type), optional, intent(inout) :: sav
597 
598  integer :: k
599  real(kind=default_precision) :: actual_gal, dtm_x2
600 
601  dtm_x2 = 2.0_default_precision * dtm
602 
603  actual_gal = merge(gal, real(0.0_default_precision, kind=default_precision), flow_field)
604 
605  sfield%data(1,y_local_index, x_local_index)=0.0_default_precision
606 
607  do k=1,local_grid%size(z_index)
608  ! Save the Z field which is used in the Robert filter
609  if (present(sav) .and. direction .eq. centred_stepping) &
610  sav%data(k,y_local_index, x_local_index) = zfield%data(k, y_local_index, x_local_index) + actual_gal
611  if (flow_field) field%data(k, y_local_index, x_local_index) = actual_gal + field%data(k, y_local_index, x_local_index)
612  if (direction == forward_stepping) then
613  zfield%data(k, y_local_index, x_local_index) = field%data(k, y_local_index, x_local_index) + dtm * &
614  sfield%data(k, y_local_index, x_local_index)
615  else
616  zfield%data(k, y_local_index, x_local_index) = actual_gal+zfield%data(k, y_local_index, x_local_index)+dtm_x2*&
617  sfield%data(k, y_local_index, x_local_index)
618  end if
619  end do
Here is the caller graph for this function:

◆ step_single_field()

subroutine stepfields_mod::step_single_field ( integer, intent(in)  x_local_index,
integer, intent(in)  y_local_index,
integer, intent(in)  x_prev,
integer, intent(in)  y_prev,
type(prognostic_field_type), intent(inout)  field,
type(prognostic_field_type), intent(inout)  zfield,
type(prognostic_field_type), intent(inout)  sfield,
type(local_grid_type), intent(inout)  local_grid,
logical, intent(in)  flow_field,
integer, intent(in)  direction,
real(kind=default_precision), intent(in)  dtm,
real(kind=default_precision), intent(in)  gal,
real(kind=default_precision), intent(in)  c1,
real(kind=default_precision), intent(in)  c2,
logical, intent(in)  do_timesmoothing,
type(prognostic_field_type), intent(inout), optional  sav 
)
private

Steps a single specific field. This will step on the yth column of the x-2 slice and x-1 and x if this is the last slice.

Parameters
x_local_indexThe current local x index
y_local_indexThe current local y index
x_prevThe previous x index to step
y_prevThe previous y index to step
fieldThe prognostic field
zfieldZ prognostic field
sfieldSource terms for the prognostic field
local_gridDescription of the local grid
flow_fieldWhether or not this is a flow field
directionThe stepping direction (centred or forward)
dtmThe delta time per timestep
galGalilean transformation
c1Constant to use in smoothing
c2Constant to use in smoothing
do_timesmoothingWhether timesmoothing using Robert filter should be done on the field
savOptional sav field

Definition at line 466 of file stepfields.F90.

468  integer, intent(in) :: x_local_index, y_local_index, x_prev, y_prev, direction
469  real(kind=default_precision), intent(in) :: dtm, gal
470  logical, intent(in) :: flow_field, do_timesmoothing
471  type(local_grid_type), intent(inout) :: local_grid
472  type(prognostic_field_type), intent(inout) :: field, zfield, sfield
473  real(kind=default_precision), intent(in) :: c1, c2
474  type(prognostic_field_type), optional, intent(inout) :: sav
475 
476  if (x_prev .ge. local_grid%local_domain_start_index(x_index)) then
477  if (present(sav)) then
478  call step_column_in_slice(y_local_index, x_prev, y_prev, field, zfield, sfield, local_grid, &
479  flow_field, direction, dtm, gal, c1, c2, do_timesmoothing, sav)
480  else
481  call step_column_in_slice(y_local_index, x_prev, y_prev, field, zfield, sfield, local_grid, &
482  flow_field, direction, dtm, gal, c1, c2, do_timesmoothing)
483  end if
484  end if
485 
486  if (x_local_index == local_grid%local_domain_end_index(x_index)) then
487  ! If this is the last slice then process x-1 (if applicable) and x too
488  if (x_local_index .gt. 1) then
489  if (present(sav)) then
490  call step_column_in_slice(y_local_index, x_local_index-1, y_prev, field, zfield, sfield, local_grid, &
491  flow_field, direction, dtm, gal, c1, c2, do_timesmoothing, sav)
492  else
493  call step_column_in_slice(y_local_index, x_local_index-1, y_prev, field, zfield, sfield, local_grid, &
494  flow_field, direction, dtm, gal, c1, c2, do_timesmoothing)
495  end if
496  end if
497  if (present(sav)) then
498  call step_column_in_slice(y_local_index, x_local_index, y_prev, field, zfield, sfield, local_grid, &
499  flow_field, direction, dtm, gal, c1, c2, do_timesmoothing, sav)
500  else
501  call step_column_in_slice(y_local_index, x_local_index, y_prev, field, zfield, sfield, local_grid, &
502  flow_field, direction, dtm, gal, c1, c2, do_timesmoothing)
503  end if
504  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ stepfields_get_descriptor()

type(component_descriptor_type) function, public stepfields_mod::stepfields_get_descriptor

Provides the descriptor back to the caller and is used in component registration.

Returns
The KidReader component descriptor

Definition at line 56 of file stepfields.F90.

57  stepfields_get_descriptor%name="stepfields"
58  stepfields_get_descriptor%version=0.1
59  stepfields_get_descriptor%initialisation=>initialisation_callback
60  stepfields_get_descriptor%timestep=>timestep_callback
61  stepfields_get_descriptor%finalisation=>finalisation_callback
62 
63  stepfields_get_descriptor%field_value_retrieval=>field_value_retrieval_callback
64  stepfields_get_descriptor%field_information_retrieval=>field_information_retrieval_callback
65  allocate(stepfields_get_descriptor%published_fields(8+8))
66 
67  stepfields_get_descriptor%published_fields(1)="tend_th_total_3d_local"
68  stepfields_get_descriptor%published_fields(2)="tend_qv_total_3d_local"
69  stepfields_get_descriptor%published_fields(3)="tend_ql_total_3d_local"
70  stepfields_get_descriptor%published_fields(4)="tend_qi_total_3d_local"
71  stepfields_get_descriptor%published_fields(5)="tend_qr_total_3d_local"
72  stepfields_get_descriptor%published_fields(6)="tend_qs_total_3d_local"
73  stepfields_get_descriptor%published_fields(7)="tend_qg_total_3d_local"
74  stepfields_get_descriptor%published_fields(8)="tend_tabs_total_3d_local"
75 
76  stepfields_get_descriptor%published_fields(8+1)="tend_th_total_profile_total_local"
77  stepfields_get_descriptor%published_fields(8+2)="tend_qv_total_profile_total_local"
78  stepfields_get_descriptor%published_fields(8+3)="tend_ql_total_profile_total_local"
79  stepfields_get_descriptor%published_fields(8+4)="tend_qi_total_profile_total_local"
80  stepfields_get_descriptor%published_fields(8+5)="tend_qr_total_profile_total_local"
81  stepfields_get_descriptor%published_fields(8+6)="tend_qs_total_profile_total_local"
82  stepfields_get_descriptor%published_fields(8+7)="tend_qg_total_profile_total_local"
83  stepfields_get_descriptor%published_fields(8+8)="tend_tabs_total_profile_total_local"
84 
Here is the call graph for this function:

◆ timestep_callback()

subroutine stepfields_mod::timestep_callback ( type(model_state_type), intent(inout), target  current_state)
private

Called at each timestep and will perform swapping and smoothing as required.

Parameters
current_stateThe current model state_mod

Definition at line 236 of file stepfields.F90.

237  type(model_state_type), target, intent(inout) :: current_state
238 
239  integer :: iq
240  integer :: current_x_index, current_y_index, target_x_index, target_y_index
241 
242  current_x_index=current_state%column_local_x
243  current_y_index=current_state%column_local_y
244  target_y_index=current_y_index-current_state%local_grid%halo_size(y_index)
245  target_x_index=current_x_index-current_state%local_grid%halo_size(x_index)
246 
247 
248  if (cfl_is_enabled .and. current_state%first_timestep_column) then
249  if (mod(current_state%timestep, current_state%cfl_frequency) == 1 .or. &
250  current_state%timestep-current_state%start_timestep .le. current_state%cfl_frequency) then
251  determine_flow_minmax=.true.
252  call reset_local_minmax_values(current_state)
253  else
254  determine_flow_minmax=.false.
255  end if
256  end if
257 
258  if (.not. current_state%halo_column) then
259  if (determine_flow_minmax .and. cfl_is_enabled) &
260  call determine_local_flow_minmax(current_state, current_state%column_local_y, current_state%column_local_x)
261  call step_all_fields(current_state)
262  end if
263 
264  ! Remove negative rounding errors
265  if (l_nonconservative_positive_q)then
266  do iq=1,current_state%number_q_fields
267  if (current_state%first_timestep_column)then
268  resetq_min(iq)=minval(current_state%zq(iq)%data(:,current_state%column_local_y, current_state%column_local_x))
269  else
270  resetq_min(iq)=min(resetq_min(iq),&
271  minval(current_state%zq(iq)%data(:,current_state%column_local_y, current_state%column_local_x)))
272  end if
273  call remove_negative_rounding_errors_for_single_field(current_state%column_local_x, current_state%column_local_y, &
274  current_state%column_local_x-2, current_state%column_local_y-1, current_state%zq(iq), current_state%local_grid)
275  end do
276  end if
277 
278 
279  ! Zero profile tendency totals on first instance in the sum
280  if (current_state%first_timestep_column) then
281  if (l_tend_pr_tot_th) then
282  tend_pr_tot_th(:)=0.0_default_precision
283  endif
284  if (l_tend_pr_tot_qv) then
285  tend_pr_tot_qv(:)=0.0_default_precision
286  endif
287  if (l_tend_pr_tot_ql) then
288  tend_pr_tot_ql(:)=0.0_default_precision
289  endif
290  if (l_tend_pr_tot_qi) then
291  tend_pr_tot_qi(:)=0.0_default_precision
292  endif
293  if (l_tend_pr_tot_qr) then
294  tend_pr_tot_qr(:)=0.0_default_precision
295  endif
296  if (l_tend_pr_tot_qs) then
297  tend_pr_tot_qs(:)=0.0_default_precision
298  endif
299  if (l_tend_pr_tot_qg) then
300  tend_pr_tot_qg(:)=0.0_default_precision
301  endif
302  if (l_tend_pr_tot_tabs) then
303  tend_pr_tot_tabs(:)=0.0_default_precision
304  endif
305  endif ! zero totals
306 
307  if (mod(current_state%timestep, diagnostic_generation_frequency) == 0 .and. .not. current_state%halo_column) then
308  call compute_component_tendencies(current_state, current_x_index, current_y_index, target_x_index, target_y_index)
309  end if
310 
Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ cfl_is_enabled

logical stepfields_mod::cfl_is_enabled
private

Definition at line 22 of file stepfields.F90.

◆ determine_flow_minmax

logical stepfields_mod::determine_flow_minmax =.false.
private

Definition at line 22 of file stepfields.F90.

22  logical :: determine_flow_minmax=.false., cfl_is_enabled

◆ diagnostic_generation_frequency

integer stepfields_mod::diagnostic_generation_frequency
private

Definition at line 47 of file stepfields.F90.

47  integer :: diagnostic_generation_frequency

◆ iqg

integer stepfields_mod::iqg =0
private

Definition at line 46 of file stepfields.F90.

◆ iqi

integer stepfields_mod::iqi =0
private

Definition at line 46 of file stepfields.F90.

◆ iql

integer stepfields_mod::iql =0
private

Definition at line 46 of file stepfields.F90.

◆ iqr

integer stepfields_mod::iqr =0
private

Definition at line 46 of file stepfields.F90.

◆ iqs

integer stepfields_mod::iqs =0
private

Definition at line 46 of file stepfields.F90.

◆ iqv

integer stepfields_mod::iqv =0
private

Definition at line 46 of file stepfields.F90.

46  integer :: iqv=0, iql=0, iqr=0, iqi=0, iqs=0, iqg=0

◆ l_nonconservative_positive_q

logical stepfields_mod::l_nonconservative_positive_q =.true.
private

Definition at line 26 of file stepfields.F90.

26  logical :: l_nonconservative_positive_q=.true.

◆ l_tend_3d_qg

logical stepfields_mod::l_tend_3d_qg
private

Definition at line 34 of file stepfields.F90.

◆ l_tend_3d_qi

logical stepfields_mod::l_tend_3d_qi
private

Definition at line 34 of file stepfields.F90.

◆ l_tend_3d_ql

logical stepfields_mod::l_tend_3d_ql
private

Definition at line 34 of file stepfields.F90.

◆ l_tend_3d_qr

logical stepfields_mod::l_tend_3d_qr
private

Definition at line 34 of file stepfields.F90.

◆ l_tend_3d_qs

logical stepfields_mod::l_tend_3d_qs
private

Definition at line 34 of file stepfields.F90.

◆ l_tend_3d_qv

logical stepfields_mod::l_tend_3d_qv
private

Definition at line 34 of file stepfields.F90.

◆ l_tend_3d_tabs

logical stepfields_mod::l_tend_3d_tabs
private

Definition at line 34 of file stepfields.F90.

◆ l_tend_3d_th

logical stepfields_mod::l_tend_3d_th
private

Definition at line 34 of file stepfields.F90.

34  logical :: l_tend_3d_th,l_tend_3d_qv, &
35  l_tend_3d_ql,l_tend_3d_qi,l_tend_3d_qr,l_tend_3d_qs,l_tend_3d_qg, &
36  l_tend_3d_tabs

◆ l_tend_pr_tot_qg

logical stepfields_mod::l_tend_pr_tot_qg
private

Definition at line 42 of file stepfields.F90.

◆ l_tend_pr_tot_qi

logical stepfields_mod::l_tend_pr_tot_qi
private

Definition at line 42 of file stepfields.F90.

◆ l_tend_pr_tot_ql

logical stepfields_mod::l_tend_pr_tot_ql
private

Definition at line 42 of file stepfields.F90.

◆ l_tend_pr_tot_qr

logical stepfields_mod::l_tend_pr_tot_qr
private

Definition at line 42 of file stepfields.F90.

◆ l_tend_pr_tot_qs

logical stepfields_mod::l_tend_pr_tot_qs
private

Definition at line 42 of file stepfields.F90.

◆ l_tend_pr_tot_qv

logical stepfields_mod::l_tend_pr_tot_qv
private

Definition at line 42 of file stepfields.F90.

◆ l_tend_pr_tot_tabs

logical stepfields_mod::l_tend_pr_tot_tabs
private

Definition at line 42 of file stepfields.F90.

◆ l_tend_pr_tot_th

logical stepfields_mod::l_tend_pr_tot_th
private

Definition at line 42 of file stepfields.F90.

42  logical :: l_tend_pr_tot_th,l_tend_pr_tot_qv, &
43  l_tend_pr_tot_ql,l_tend_pr_tot_qi,l_tend_pr_tot_qr,l_tend_pr_tot_qs,l_tend_pr_tot_qg, &
44  l_tend_pr_tot_tabs

◆ resetq_min

real(kind=default_precision), dimension(:), allocatable stepfields_mod::resetq_min
private

Definition at line 25 of file stepfields.F90.

25  real(kind=default_precision), allocatable :: resetq_min(:)

◆ tend_3d_qg

real(kind=default_precision), dimension(:,:,:), allocatable stepfields_mod::tend_3d_qg
private

Definition at line 30 of file stepfields.F90.

◆ tend_3d_qi

real(kind=default_precision), dimension(:,:,:), allocatable stepfields_mod::tend_3d_qi
private

Definition at line 30 of file stepfields.F90.

◆ tend_3d_ql

real(kind=default_precision), dimension(:,:,:), allocatable stepfields_mod::tend_3d_ql
private

Definition at line 30 of file stepfields.F90.

◆ tend_3d_qr

real(kind=default_precision), dimension(:,:,:), allocatable stepfields_mod::tend_3d_qr
private

Definition at line 30 of file stepfields.F90.

◆ tend_3d_qs

real(kind=default_precision), dimension(:,:,:), allocatable stepfields_mod::tend_3d_qs
private

Definition at line 30 of file stepfields.F90.

◆ tend_3d_qv

real(kind=default_precision), dimension(:,:,:), allocatable stepfields_mod::tend_3d_qv
private

Definition at line 30 of file stepfields.F90.

◆ tend_3d_tabs

real(kind=default_precision), dimension(:,:,:), allocatable stepfields_mod::tend_3d_tabs
private

Definition at line 30 of file stepfields.F90.

◆ tend_3d_th

real(kind=default_precision), dimension(:,:,:), allocatable stepfields_mod::tend_3d_th
private

Definition at line 30 of file stepfields.F90.

30  real(kind=default_precision), dimension(:,:,:), allocatable :: &
31  tend_3d_th,tend_3d_qv, &
32  tend_3d_ql,tend_3d_qi,tend_3d_qr,tend_3d_qs,tend_3d_qg, &
33  tend_3d_tabs

◆ tend_pr_tot_qg

real(kind=default_precision), dimension(:), allocatable stepfields_mod::tend_pr_tot_qg
private

Definition at line 38 of file stepfields.F90.

◆ tend_pr_tot_qi

real(kind=default_precision), dimension(:), allocatable stepfields_mod::tend_pr_tot_qi
private

Definition at line 38 of file stepfields.F90.

◆ tend_pr_tot_ql

real(kind=default_precision), dimension(:), allocatable stepfields_mod::tend_pr_tot_ql
private

Definition at line 38 of file stepfields.F90.

◆ tend_pr_tot_qr

real(kind=default_precision), dimension(:), allocatable stepfields_mod::tend_pr_tot_qr
private

Definition at line 38 of file stepfields.F90.

◆ tend_pr_tot_qs

real(kind=default_precision), dimension(:), allocatable stepfields_mod::tend_pr_tot_qs
private

Definition at line 38 of file stepfields.F90.

◆ tend_pr_tot_qv

real(kind=default_precision), dimension(:), allocatable stepfields_mod::tend_pr_tot_qv
private

Definition at line 38 of file stepfields.F90.

◆ tend_pr_tot_tabs

real(kind=default_precision), dimension(:), allocatable stepfields_mod::tend_pr_tot_tabs
private

Definition at line 38 of file stepfields.F90.

◆ tend_pr_tot_th

real(kind=default_precision), dimension(:), allocatable stepfields_mod::tend_pr_tot_th
private

Definition at line 38 of file stepfields.F90.

38  real(kind=default_precision), dimension(:), allocatable :: &
39  tend_pr_tot_th,tend_pr_tot_qv, &
40  tend_pr_tot_ql,tend_pr_tot_qi,tend_pr_tot_qr,tend_pr_tot_qs,tend_pr_tot_qg, &
41  tend_pr_tot_tabs
grids_mod::x_index
integer, parameter, public x_index
Definition: grids.F90:14
optionsdatabase_mod::options_get_integer
integer function, public options_get_integer(options_database, key, index)
Retrieves an integer value from the database that matches the provided key.
Definition: optionsdatabase.F90:217
grids_mod::y_index
integer, parameter, public y_index
Definition: grids.F90:14
grids_mod::z_index
integer, parameter, public z_index
Grid index parameters.
Definition: grids.F90:14
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