MONC
Functions/Subroutines | Variables
cfltest_mod Module Reference

This contains the CFL test. It will perform the local advective CFL and Galilean transfromation calculations, compute the global values of these, then check the cfl criterion and determine an absolute new dtm. Depending upon the maximum and increment values this is then smoothed into a new dtm value. The value of dtm is physically set to this new value at the start of the next timestep. More...

Functions/Subroutines

type(component_descriptor_type) function, public cfltest_get_descriptor ()
 Provides the descriptor back to the caller and is used in component registration. More...
 
subroutine initialisation_callback (current_state)
 Called at initialisation, will read in configuration and use either configured or default values. More...
 
subroutine timestep_callback (current_state)
 Called at each timestep, this will only do the CFL computation every nncfl timesteps (or every timestep up to nncfl) but will ratchet up to the absolute (target) dtm as needed. More...
 
subroutine update_dtm_based_on_absolute (current_state, cfl_number)
 Updates the (new) dtm value, which is actioned after time step completion, based upon the absolute value. This is incremented towards the absolute if that is too large a step, and even if the absolute value has not been updated in this timestep, this ratcheting will still occur if needed. More...
 
subroutine perform_cfl_and_galilean_transformation_calculation (current_state)
 Performs the CFL and Galilean transformation calculations. First locally and then will determine the global value of each calculation. If U, V or W are not active then these are set to 0 and the calculation use these values. More...
 
subroutine get_global_values (local_zumin, local_zumax, local_zvmin, local_zvmax, local_cvel_z, local_cvis, global_zumin, global_zumax, global_zvmin, global_zvmax, global_cvel_z, global_cvis, parallel_state)
 Gets the global reduction values based upon the local contributions of CFL and Galilean transformations provided. More...
 

Variables

real(kind=default_precision) tollerance
 
real(kind=default_precision) cvismax
 
real(kind=default_precision) cvelmax
 
real(kind=default_precision) dtmmax
 
real(kind=default_precision) dtmmin
 
real(kind=default_precision) rincmax
 
logical l_monitor_cfl
 

Detailed Description

This contains the CFL test. It will perform the local advective CFL and Galilean transfromation calculations, compute the global values of these, then check the cfl criterion and determine an absolute new dtm. Depending upon the maximum and increment values this is then smoothed into a new dtm value. The value of dtm is physically set to this new value at the start of the next timestep.

Function/Subroutine Documentation

◆ cfltest_get_descriptor()

type(component_descriptor_type) function, public cfltest_mod::cfltest_get_descriptor

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

Returns
The termination check component descriptor

Definition at line 31 of file cfltest.F90.

32  cfltest_get_descriptor%name="cfltest"
33  cfltest_get_descriptor%version=0.1
34  cfltest_get_descriptor%initialisation=>initialisation_callback
35  cfltest_get_descriptor%timestep=>timestep_callback
Here is the call graph for this function:

◆ get_global_values()

subroutine cfltest_mod::get_global_values ( real(kind=default_precision), intent(in)  local_zumin,
real(kind=default_precision), intent(in)  local_zumax,
real(kind=default_precision), intent(in)  local_zvmin,
real(kind=default_precision), intent(in)  local_zvmax,
real(kind=default_precision), intent(in)  local_cvel_z,
real(kind=default_precision), intent(in)  local_cvis,
real(kind=default_precision), intent(out)  global_zumin,
real(kind=default_precision), intent(out)  global_zumax,
real(kind=default_precision), intent(out)  global_zvmin,
real(kind=default_precision), intent(out)  global_zvmax,
real(kind=default_precision), intent(out)  global_cvel_z,
real(kind=default_precision), intent(out)  global_cvis,
type(parallel_state_type), intent(inout)  parallel_state 
)
private

Gets the global reduction values based upon the local contributions of CFL and Galilean transformations provided.

Parameters
local_zuminLocal contribution to ZU minimum
local_zumaxLocal contribution to ZU maximum
local_zvminLocal contribution to ZV minimum
local_zvmaxLocal contribution to ZV maximum
local_cvel_zLocal contribution to component of Courant number in z
local_cvisLocal contribution to viscous Courant number
global_zuminGlobally reduced value of ZU minimum
global_zumaxGlobally reduced value of ZU maximum
global_zvminGlobally reduced value of ZV minimum
global_zvmaxGlobally reduced value of ZV maximum
global_cvel_zGlobally reduced value of contribution to component of Courant number in z
global_cvisGlobally reduced value of contribution to viscous Courant number
parallel_stateThe parallel state

Definition at line 200 of file cfltest.F90.

202  type(parallel_state_type), intent(inout) :: parallel_state
203  real(kind=default_precision), intent(in) :: local_zumin, local_zumax, local_zvmin, local_zvmax, local_cvel_z, local_cvis
204  real(kind=default_precision), intent(out) :: global_zumin, global_zumax, global_zvmin, global_zvmax, global_cvel_z, global_cvis
205 
206  integer :: ierr
207 
208  call mpi_allreduce(local_zumax, global_zumax, 1, precision_type, mpi_max, parallel_state%monc_communicator, ierr)
209  call mpi_allreduce(local_zvmax, global_zvmax, 1, precision_type, mpi_max, parallel_state%monc_communicator, ierr)
210  call mpi_allreduce(local_cvel_z, global_cvel_z, 1, precision_type, mpi_max, parallel_state%monc_communicator, ierr)
211  call mpi_allreduce(local_cvis, global_cvis, 1, precision_type, mpi_max, parallel_state%monc_communicator, ierr)
212  call mpi_allreduce(local_zumin, global_zumin, 1, precision_type, mpi_min, parallel_state%monc_communicator, ierr)
213  call mpi_allreduce(local_zvmin, global_zvmin, 1, precision_type, mpi_min, parallel_state%monc_communicator, ierr)
Here is the caller graph for this function:

◆ initialisation_callback()

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

Called at initialisation, will read in configuration and use either configured or default values.

Parameters
current_stateThe current model state

Definition at line 40 of file cfltest.F90.

41  type(model_state_type), intent(inout), target :: current_state
42 
43  current_state%cfl_frequency=options_get_integer(current_state%options_database, "cfl_frequency")
44  tollerance=options_get_real(current_state%options_database, "cfl_tollerance")
45  cvismax=options_get_real(current_state%options_database, "cfl_cvismax")
46  cvelmax=options_get_real(current_state%options_database, "cfl_cvelmax")
47  dtmmax=options_get_real(current_state%options_database, "cfl_dtmmax")
48  dtmmin=options_get_real(current_state%options_database, "cfl_dtmmin")
49  rincmax=options_get_real(current_state%options_database, "cfl_rincmax")
50 
51  l_monitor_cfl = options_get_logical(current_state%options_database,"cfl_monitor")
52 
53  allocate(current_state%abswmax(current_state%local_grid%local_domain_end_index(z_index)))
Here is the caller graph for this function:

◆ perform_cfl_and_galilean_transformation_calculation()

subroutine cfltest_mod::perform_cfl_and_galilean_transformation_calculation ( type(model_state_type), intent(inout), target  current_state)
private

Performs the CFL and Galilean transformation calculations. First locally and then will determine the global value of each calculation. If U, V or W are not active then these are set to 0 and the calculation use these values.

Parameters
current_stateThe current model state

Definition at line 138 of file cfltest.F90.

139  type(model_state_type), intent(inout), target :: current_state
140 
141  integer :: k
142  real(kind=default_precision) :: global_zumin, global_zumax, global_zvmin, &
143  global_zvmax, global_cvel_z, global_cvis
144 
145 #ifdef U_ACTIVE
146  current_state%local_zumin=current_state%local_zumin+current_state%ugal ! _undo Gal-trfm
147  current_state%local_zumax=current_state%local_zumax+current_state%ugal
148 #else
149  current_state%local_zumin=0.0_default_precision
150  current_state%local_zumax=0.0_default_precision
151 #endif
152 #ifdef V_ACTIVE
153  current_state%local_zvmin=current_state%local_zvmin+current_state%vgal ! _undo Gal-trfm
154  current_state%local_zvmax=current_state%local_zvmax+current_state%vgal
155 #else
156  current_state%local_zvmin=0.0_default_precision
157  current_state%local_zvmax=0.0_default_precision
158 #endif
159 #ifdef W_ACTIVE
160  current_state%local_cvel_z=current_state%cvel_z
161  do k=2,current_state%local_grid%local_domain_end_index(z_index)-1
162  ! CVELZ will be multiplied by DTM in TESTCFL
163  current_state%local_cvel_z=max(current_state%local_cvel_z, &
164  current_state%abswmax(k)*current_state%global_grid%configuration%vertical%rdzn(k+1))
165  end do
166 #else
167  current_state%local_cvel_z=0.0_default_precision
168 #endif
169  call get_global_values(current_state%local_zumin, current_state%local_zumax, current_state%local_zvmin, &
170  current_state%local_zvmax, current_state%local_cvel_z, current_state%cvis, &
171  global_zumin, global_zumax, global_zvmin, global_zvmax, global_cvel_z, global_cvis, current_state%parallel)
172 
173  if (current_state%galilean_transformation) then
174  if (.not.current_state%fix_ugal)current_state%ugal=0.5_default_precision*(global_zumin+global_zumax)
175  if (.not.current_state%fix_vgal)current_state%vgal=0.5_default_precision*(global_zvmin+global_zvmax)
176  else
177  current_state%ugal=0.0_default_precision
178  current_state%vgal=0.0_default_precision
179  end if
180  current_state%cvel_z=global_cvel_z
181  current_state%cvel_x=max(abs(global_zumax-current_state%ugal), abs(global_zumin-current_state%ugal))
182  current_state%cvel_y=max(abs(global_zvmax-current_state%vgal), abs(global_zvmin-current_state%vgal))
183  current_state%cvis=global_cvis
Here is the call graph for this function:
Here is the caller graph for this function:

◆ timestep_callback()

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

Called at each timestep, this will only do the CFL computation every nncfl timesteps (or every timestep up to nncfl) but will ratchet up to the absolute (target) dtm as needed.

Parameters
current_stateThe current model state

Definition at line 59 of file cfltest.F90.

60  type(model_state_type), intent(inout), target :: current_state
61 
62  real(kind=default_precision) :: cfl_number
63 
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  current_state%cvel=0.0_default_precision
67  current_state%cvel_x=0.0_default_precision
68  current_state%cvel_y=0.0_default_precision
69  current_state%cvel_z=0.0_default_precision
70 
71  call perform_cfl_and_galilean_transformation_calculation(current_state)
72 
73  current_state%cvel=(current_state%cvel_x*current_state%global_grid%configuration%horizontal%cx+current_state%cvel_y*&
74  current_state%global_grid%configuration%horizontal%cy+current_state%cvel_z)*current_state%dtm
75  current_state%cvis=current_state%cvis*(current_state%dtm * 4)
76 
77  cfl_number=current_state%cvis/cvismax+current_state%cvel/cvelmax
78 
79  current_state%absolute_new_dtm=current_state%dtm
80  current_state%update_dtm=.false.
81  if (cfl_number .gt. 0.0_default_precision) then
82  if (cfl_number .lt. (1.0_default_precision-tollerance) .or. cfl_number .gt. (1.0_default_precision+tollerance)) then
83  current_state%absolute_new_dtm=current_state%dtm/cfl_number
84  end if
85  end if
86  end if
87  call update_dtm_based_on_absolute(current_state, cfl_number)
88  current_state%cvis=0.0_default_precision
Here is the call graph for this function:
Here is the caller graph for this function:

◆ update_dtm_based_on_absolute()

subroutine cfltest_mod::update_dtm_based_on_absolute ( type(model_state_type), intent(inout), target  current_state,
real(kind=default_precision), intent(in)  cfl_number 
)
private

Updates the (new) dtm value, which is actioned after time step completion, based upon the absolute value. This is incremented towards the absolute if that is too large a step, and even if the absolute value has not been updated in this timestep, this ratcheting will still occur if needed.

Parameters
current_stateThe current model state

Definition at line 95 of file cfltest.F90.

96  type(model_state_type), intent(inout), target :: current_state
97  real(kind=default_precision), intent(in) :: cfl_number
98 
99  if (current_state%dtm .ne. current_state%absolute_new_dtm .and. &
100  (current_state%dtm .ne. dtmmax .or. current_state%absolute_new_dtm .lt. dtmmax)) then
101 
102  current_state%update_dtm=.true.
103 
104  current_state%dtm_new=min(current_state%dtm*(1.0_default_precision+rincmax), current_state%absolute_new_dtm, dtmmax)
105 
106  !! --- Diagnostic Writing -----------------
107  if (current_state%parallel%my_rank==0) then
108  if (log_get_logging_level() .eq. log_debug) then
109  call log_log(log_debug, "dtm changed from "//trim(conv_to_string(current_state%dtm, 5))//" to "//&
110  trim(conv_to_string(current_state%dtm_new, 5)))
111  end if
112  if (current_state%dtm_new .lt. dtmmin) then
113  call log_log(log_error, "Timestep too small, dtmnew="//trim(conv_to_string(current_state%dtm_new, 5))//&
114  " dtmmin="//trim(conv_to_string(dtmmin, 5)))
115  end if
116  if (l_monitor_cfl) then
117  call log_log(log_info, " --- CFL Monitoring Information --- ")
118  call log_log(log_info, "dtm changed from "//trim(conv_to_string(current_state%dtm, 5))//" to "//&
119  trim(conv_to_string(current_state%dtm_new, 5)))
120  if (cfl_number .gt. 0.0) then
121  call log_log(log_info, "cfl_number : "//trim(conv_to_string(cfl_number))//" (change divisor)")
122  call log_log(log_info, "cvis : "//trim(conv_to_string(current_state%cvis)) )
123  call log_log(log_info, "cvel : "//trim(conv_to_string(current_state%cvel)) )
124  else
125  call log_log(log_info, "dtm change due to ratcheting only. Target dtm unchanged.")
126  end if
127  call log_log(log_info, "target dtm : "//trim(conv_to_string(current_state%absolute_new_dtm)) )
128  call log_newline()
129 
130  end if ! l_monitor_cfl
131  end if ! Diagnostic Writing
132  end if
Here is the caller graph for this function:

Variable Documentation

◆ cvelmax

real(kind=default_precision) cfltest_mod::cvelmax
private

Definition at line 23 of file cfltest.F90.

◆ cvismax

real(kind=default_precision) cfltest_mod::cvismax
private

Definition at line 23 of file cfltest.F90.

◆ dtmmax

real(kind=default_precision) cfltest_mod::dtmmax
private

Definition at line 23 of file cfltest.F90.

◆ dtmmin

real(kind=default_precision) cfltest_mod::dtmmin
private

Definition at line 23 of file cfltest.F90.

◆ l_monitor_cfl

logical cfltest_mod::l_monitor_cfl
private

Definition at line 24 of file cfltest.F90.

24  logical l_monitor_cfl

◆ rincmax

real(kind=default_precision) cfltest_mod::rincmax
private

Definition at line 23 of file cfltest.F90.

◆ tollerance

real(kind=default_precision) cfltest_mod::tollerance
private

Definition at line 23 of file cfltest.F90.

23  real(kind=default_precision) :: tollerance, cvismax, cvelmax, dtmmax, dtmmin, rincmax
logging_mod::log_error
integer, parameter, public log_error
Only log ERROR messages.
Definition: logging.F90:11
logging_mod::log_info
integer, parameter, public log_info
Log INFO, WARNING and ERROR messages.
Definition: logging.F90:13
logging_mod::log_log
subroutine, public log_log(level, message, str)
Logs a message at the specified level. If the level is above the current level then the message is ig...
Definition: logging.F90:75
logging_mod::log_get_logging_level
integer function, public log_get_logging_level()
Retrieves the current logging level.
Definition: logging.F90:122
logging_mod::log_debug
integer, parameter, public log_debug
Log DEBUG, INFO, WARNING and ERROR messages.
Definition: logging.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