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...
 

Variables

logical determine_flow_minmax =.false.
 
logical cfl_is_enabled
 
real(kind=default_precision), dimension(:), allocatable resetq_min
 
logical l_nonconservative_positive_q =.true.
 

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

◆ 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 146 of file stepfields.F90.

146  type(model_state_type), target, intent(inout) :: current_state
147  integer, intent(in) :: local_x, local_y
148 
149  integer :: k
150 
151  do k=2, current_state%local_grid%local_domain_end_index(z_index)
152 #ifdef U_ACTIVE
153  current_state%local_zumax = max(current_state%local_zumax, current_state%zu%data(k,local_y,local_x))
154  current_state%local_zumin = min(current_state%local_zumin, current_state%zu%data(k,local_y,local_x))
155 #endif
156 #ifdef V_ACTIVE
157  current_state%local_zvmax = max(current_state%local_zvmax, current_state%zv%data(k,local_y,local_x))
158  current_state%local_zvmin = min(current_state%local_zvmin, current_state%zv%data(k,local_y,local_x))
159 #endif
160 #ifdef W_ACTIVE
161  if (k .lt. current_state%local_grid%local_domain_end_index(z_index)) then
162  current_state%abswmax(k) = max(current_state%abswmax(k), abs(current_state%zw%data(k,local_y,local_x)))
163  end if
164 #endif
165  end do
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 51 of file stepfields.F90.

51  type(model_state_type), target, intent(inout) :: current_state
52 
53  deallocate(resetq_min)
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 40 of file stepfields.F90.

40  type(model_state_type), target, intent(inout) :: current_state
41 
42  allocate(resetq_min(current_state%number_q_fields))
43  cfl_is_enabled=is_component_enabled(current_state%options_database, "cfltest")
44  if (cfl_is_enabled) call reset_local_minmax_values(current_state)
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 296 of file stepfields.F90.

296  type(prognostic_field_type), intent(inout) :: field, zfield
297  type(local_grid_type), intent(inout) :: local_grid
298  integer, intent(in) :: x_index, y_index
299  real(kind=DEFAULT_PRECISION), intent(in) :: c1, c2
300 
301  integer :: k
302 
303  do k=1,local_grid%size(z_index)
304  field%data(k, y_index, x_index)=c1*field%data(k, y_index, x_index)+c2*zfield%data(k, y_index, x_index)
305  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 190 of file stepfields.F90.

190  integer, intent(in) :: x_local_index, y_local_index, x_prev, y_prev
191  type(local_grid_type), intent(inout) :: local_grid
192  type(prognostic_field_type), intent(inout) :: field
193 
194  if (x_prev .ge. local_grid%local_domain_start_index(x_index)) then
195  call remove_negative_rounding_errors_in_slice(y_local_index, x_prev, y_prev, field, local_grid)
196  end if
197  if (x_local_index == local_grid%local_domain_end_index(x_index)) then
198  if (x_local_index .gt. 1) then
199  call remove_negative_rounding_errors_in_slice(y_local_index, x_local_index-1, y_prev, field, local_grid)
200  end if
201  call remove_negative_rounding_errors_in_slice(y_local_index, x_local_index, y_prev, field, local_grid)
202  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 213 of file stepfields.F90.

213  integer, intent(in) :: y_local_index, x_prev, y_prev
214  type(local_grid_type), intent(inout) :: local_grid
215  type(prognostic_field_type), intent(inout) :: field
216 
217  if (y_prev .ge. local_grid%local_domain_start_index(y_index)) then
218  where (field%data(:, y_prev, x_prev) < 0.0_default_precision)
219  field%data(:, y_prev, x_prev)=0.0_default_precision
220  end where
221  end if
222  if (y_local_index == local_grid%local_domain_end_index(y_index)) then
223  where (field%data(:, y_local_index, x_prev) < 0.0_default_precision)
224  field%data(:, y_local_index, x_prev)=0.0_default_precision
225  end where
226  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 171 of file stepfields.F90.

171  type(model_state_type), intent(inout), target :: current_state
172 
173  ! Reset the local values for the next timestep
174  current_state%local_zumin=rlargep
175  current_state%local_zumax=-rlargep
176  current_state%local_zvmin=rlargep
177  current_state%local_zvmax=-rlargep
178  current_state%abswmax=-rlargep
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 98 of file stepfields.F90.

98  type(model_state_type), target, intent(inout) :: current_state
99 
100  integer :: x_prev, y_prev, i, k
101  real(kind=DEFAULT_PRECISION) :: c1, c2
102 
103  x_prev = current_state%column_local_x-2
104  y_prev = current_state%column_local_y-1
105 
106  c1 = 1.0_default_precision - 2.0_default_precision*current_state%tsmth
107  c2 = current_state%tsmth
108 
109 #ifdef U_ACTIVE
110  call step_single_field(current_state%column_local_x, current_state%column_local_y, &
111  x_prev, y_prev, current_state%u, current_state%zu, current_state%su, current_state%local_grid, .true., &
112  current_state%field_stepping, current_state%dtm, current_state%ugal, c1, c2, .false., current_state%savu)
113 #endif
114 #ifdef V_ACTIVE
115  call step_single_field(current_state%column_local_x, current_state%column_local_y, &
116  x_prev, y_prev, current_state%v, current_state%zv, current_state%sv, current_state%local_grid, .true., &
117  current_state%field_stepping, current_state%dtm, current_state%vgal, c1, c2, .false., current_state%savv)
118 #endif
119 #ifdef W_ACTIVE
120  call step_single_field(current_state%column_local_x, current_state%column_local_y, &
121  x_prev, y_prev, current_state%w, current_state%zw, current_state%sw, current_state%local_grid, .false., &
122  current_state%field_stepping, current_state%dtm, real(0., kind=DEFAULT_PRECISION), c1, c2, .false., current_state%savw)
123 #endif
124  if (current_state%th%active) then
125  call step_single_field(current_state%column_local_x, current_state%column_local_y, &
126  x_prev, y_prev, current_state%th, current_state%zth, current_state%sth, current_state%local_grid, .false., &
127  current_state%field_stepping, current_state%dtm, real(0., kind=DEFAULT_PRECISION), c1, c2, &
128  current_state%field_stepping == centred_stepping)
129  endif
130  do i=1,current_state%number_q_fields
131  if (current_state%q(i)%active) then
132  call step_single_field(current_state%column_local_x, current_state%column_local_y, x_prev, y_prev, &
133  current_state%q(i), current_state%zq(i), current_state%sq(i), current_state%local_grid, .false., &
134  current_state%field_stepping, current_state%dtm, real(0., kind=DEFAULT_PRECISION), c1, c2, &
135  current_state%field_stepping == centred_stepping)
136  end if
137  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 328 of file stepfields.F90.

328  integer, intent(in) :: y_local_index, x_prev, y_prev, direction
329  real(kind=DEFAULT_PRECISION), intent(in) :: dtm, gal
330  logical, intent(in) :: flow_field, do_timesmoothing
331  type(local_grid_type), intent(inout) :: local_grid
332  type(prognostic_field_type), intent(inout) :: field, zfield, sfield
333  real(kind=DEFAULT_PRECISION), intent(in) :: c1, c2
334  type(prognostic_field_type), optional, intent(inout) :: sav
335 
336  if (y_prev .ge. local_grid%local_domain_start_index(y_index)) then
337  if (do_timesmoothing) then
338  call perform_timesmooth_for_field(field, zfield, local_grid, x_prev, y_prev, c1, c2)
339  end if
340  if (present(sav)) then
341  call step_field(x_prev, y_prev, field, zfield, sfield, local_grid, flow_field, direction, dtm, gal, sav)
342  else
343  call step_field(x_prev, y_prev, field, zfield, sfield, local_grid, flow_field, direction, dtm, gal)
344  end if
345  end if
346 
347  if (y_local_index == local_grid%local_domain_end_index(y_index)) then
348  if (do_timesmoothing) then
349  call perform_timesmooth_for_field(field, zfield, local_grid, x_prev, y_local_index, c1, c2)
350  end if
351  if (present(sav)) then
352  call step_field(x_prev, y_local_index, field, zfield, sfield, local_grid, flow_field, direction, dtm, gal, sav)
353  else
354  call step_field(x_prev, y_local_index, field, zfield, sfield, local_grid, flow_field, direction, dtm, gal)
355  end if
356  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 371 of file stepfields.F90.

371  integer, intent(in) :: x_local_index, y_local_index, direction
372  real(kind=DEFAULT_PRECISION), intent(in) :: dtm, gal
373  logical, intent(in) :: flow_field
374  type(local_grid_type), intent(inout) :: local_grid
375  type(prognostic_field_type), intent(inout) :: field, zfield, sfield
376  type(prognostic_field_type), optional, intent(inout) :: sav
377 
378  integer :: k
379  real(kind=DEFAULT_PRECISION) :: actual_gal, dtm_x2
380 
381  dtm_x2 = 2.0_default_precision * dtm
382 
383  actual_gal = merge(gal, real(0.0_DEFAULT_PRECISION, kind=DEFAULT_PRECISION), flow_field)
384 
385  sfield%data(1,y_local_index, x_local_index)=0.0_default_precision
386 
387  do k=1,local_grid%size(z_index)
388  ! Save the Z field which is used in the Robert filter
389  if (present(sav) .and. direction .eq. centred_stepping) &
390  sav%data(k,y_local_index, x_local_index) = zfield%data(k, y_local_index, x_local_index) + actual_gal
391  if (flow_field) field%data(k, y_local_index, x_local_index) = actual_gal + field%data(k, y_local_index, x_local_index)
392  if (direction == forward_stepping) then
393  zfield%data(k, y_local_index, x_local_index) = field%data(k, y_local_index, x_local_index) + dtm * &
394  sfield%data(k, y_local_index, x_local_index)
395  else
396  zfield%data(k, y_local_index, x_local_index) = actual_gal+zfield%data(k, y_local_index, x_local_index)+dtm_x2*&
397  sfield%data(k, y_local_index, x_local_index)
398  end if
399  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 248 of file stepfields.F90.

248  integer, intent(in) :: x_local_index, y_local_index, x_prev, y_prev, direction
249  real(kind=DEFAULT_PRECISION), intent(in) :: dtm, gal
250  logical, intent(in) :: flow_field, do_timesmoothing
251  type(local_grid_type), intent(inout) :: local_grid
252  type(prognostic_field_type), intent(inout) :: field, zfield, sfield
253  real(kind=DEFAULT_PRECISION), intent(in) :: c1, c2
254  type(prognostic_field_type), optional, intent(inout) :: sav
255 
256  if (x_prev .ge. local_grid%local_domain_start_index(x_index)) then
257  if (present(sav)) then
258  call step_column_in_slice(y_local_index, x_prev, y_prev, field, zfield, sfield, local_grid, &
259  flow_field, direction, dtm, gal, c1, c2, do_timesmoothing, sav)
260  else
261  call step_column_in_slice(y_local_index, x_prev, y_prev, field, zfield, sfield, local_grid, &
262  flow_field, direction, dtm, gal, c1, c2, do_timesmoothing)
263  end if
264  end if
265 
266  if (x_local_index == local_grid%local_domain_end_index(x_index)) then
267  ! If this is the last slice then process x-1 (if applicable) and x too
268  if (x_local_index .gt. 1) then
269  if (present(sav)) then
270  call step_column_in_slice(y_local_index, x_local_index-1, y_prev, field, zfield, sfield, local_grid, &
271  flow_field, direction, dtm, gal, c1, c2, do_timesmoothing, sav)
272  else
273  call step_column_in_slice(y_local_index, x_local_index-1, y_prev, field, zfield, sfield, local_grid, &
274  flow_field, direction, dtm, gal, c1, c2, do_timesmoothing)
275  end if
276  end if
277  if (present(sav)) then
278  call step_column_in_slice(y_local_index, x_local_index, y_prev, field, zfield, sfield, local_grid, &
279  flow_field, direction, dtm, gal, c1, c2, do_timesmoothing, sav)
280  else
281  call step_column_in_slice(y_local_index, x_local_index, y_prev, field, zfield, sfield, local_grid, &
282  flow_field, direction, dtm, gal, c1, c2, do_timesmoothing)
283  end if
284  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 30 of file stepfields.F90.

30  stepfields_get_descriptor%name="stepfields"
31  stepfields_get_descriptor%version=0.1
32  stepfields_get_descriptor%initialisation=>initialisation_callback
33  stepfields_get_descriptor%timestep=>timestep_callback
34  stepfields_get_descriptor%finalisation=>finalisation_callback
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 59 of file stepfields.F90.

59  type(model_state_type), target, intent(inout) :: current_state
60 
61  integer :: iq
62 
63  if (cfl_is_enabled .and. current_state%first_timestep_column) then
64  if (mod(current_state%timestep, current_state%cfl_frequency) == 1 .or. &
65  current_state%timestep-current_state%start_timestep .le. current_state%cfl_frequency) then
66  determine_flow_minmax=.true.
67  call reset_local_minmax_values(current_state)
68  else
69  determine_flow_minmax=.false.
70  end if
71  end if
72 
73  if (.not. current_state%halo_column) then
74  if (determine_flow_minmax .and. cfl_is_enabled) &
75  call determine_local_flow_minmax(current_state, current_state%column_local_y, current_state%column_local_x)
76  call step_all_fields(current_state)
77  end if
78 
79  ! Remove negative rounding errors
80  if (l_nonconservative_positive_q)then
81  do iq=1,current_state%number_q_fields
82  if (current_state%first_timestep_column)then
83  resetq_min(iq)=minval(current_state%zq(iq)%data(:,current_state%column_local_y, current_state%column_local_x))
84  else
85  resetq_min(iq)=min(resetq_min(iq),&
86  minval(current_state%zq(iq)%data(:,current_state%column_local_y, current_state%column_local_x)))
87  end if
88  call remove_negative_rounding_errors_for_single_field(current_state%column_local_x, current_state%column_local_y, &
89  current_state%column_local_x-2, current_state%column_local_y-1, current_state%zq(iq), current_state%local_grid)
90  end do
91  end if
92 
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 17 of file stepfields.F90.

◆ determine_flow_minmax

logical stepfields_mod::determine_flow_minmax =.false.
private

Definition at line 17 of file stepfields.F90.

17  logical :: determine_flow_minmax=.false., cfl_is_enabled

◆ l_nonconservative_positive_q

logical stepfields_mod::l_nonconservative_positive_q =.true.
private

Definition at line 21 of file stepfields.F90.

21  logical :: l_nonconservative_positive_q=.true.

◆ resetq_min

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

Definition at line 20 of file stepfields.F90.

20  real(kind=DEFAULT_PRECISION), allocatable :: resetq_min(:)