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

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 29 of file cfltest.F90.

29  cfltest_get_descriptor%name="cfltest"
30  cfltest_get_descriptor%version=0.1
31  cfltest_get_descriptor%initialisation=>initialisation_callback
32  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 177 of file cfltest.F90.

177  type(parallel_state_type), intent(inout) :: parallel_state
178  real(kind=DEFAULT_PRECISION), intent(in) :: local_zumin, local_zumax, local_zvmin, local_zvmax, local_cvel_z, local_cvis
179  real(kind=DEFAULT_PRECISION), intent(out) :: global_zumin, global_zumax, global_zvmin, global_zvmax, global_cvel_z, global_cvis
180 
181  integer :: ierr
182 
183  call mpi_allreduce(local_zumax, global_zumax, 1, precision_type, mpi_max, parallel_state%monc_communicator, ierr)
184  call mpi_allreduce(local_zvmax, global_zvmax, 1, precision_type, mpi_max, parallel_state%monc_communicator, ierr)
185  call mpi_allreduce(local_cvel_z, global_cvel_z, 1, precision_type, mpi_max, parallel_state%monc_communicator, ierr)
186  call mpi_allreduce(local_cvis, global_cvis, 1, precision_type, mpi_max, parallel_state%monc_communicator, ierr)
187  call mpi_allreduce(local_zumin, global_zumin, 1, precision_type, mpi_min, parallel_state%monc_communicator, ierr)
188  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 38 of file cfltest.F90.

38  type(model_state_type), intent(inout), target :: current_state
39 
40  current_state%cfl_frequency=options_get_integer(current_state%options_database, "cfl_frequency")
41  tollerance=options_get_real(current_state%options_database, "cfl_tollerance")
42  cvismax=options_get_real(current_state%options_database, "cfl_cvismax")
43  cvelmax=options_get_real(current_state%options_database, "cfl_cvelmax")
44  dtmmax=options_get_real(current_state%options_database, "cfl_dtmmax")
45  dtmmin=options_get_real(current_state%options_database, "cfl_dtmmin")
46  rincmax=options_get_real(current_state%options_database, "cfl_rincmax")
47 
48  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 114 of file cfltest.F90.

114  type(model_state_type), intent(inout), target :: current_state
115 
116  integer :: k
117  real(kind=DEFAULT_PRECISION) :: global_zumin, global_zumax, global_zvmin, &
118  global_zvmax, global_cvel_z, global_cvis
119 
120 #ifdef U_ACTIVE
121  current_state%local_zumin=current_state%local_zumin+current_state%ugal ! _undo Gal-trfm
122  current_state%local_zumax=current_state%local_zumax+current_state%ugal
123 #else
124  current_state%local_zumin=0.0_default_precision
125  current_state%local_zumax=0.0_default_precision
126 #endif
127 #ifdef V_ACTIVE
128  current_state%local_zvmin=current_state%local_zvmin+current_state%vgal ! _undo Gal-trfm
129  current_state%local_zvmax=current_state%local_zvmax+current_state%vgal
130 #else
131  current_state%local_zvmin=0.0_default_precision
132  current_state%local_zvmax=0.0_default_precision
133 #endif
134 #ifdef W_ACTIVE
135  current_state%local_cvel_z=current_state%cvel_z
136  do k=2,current_state%local_grid%local_domain_end_index(z_index)-1
137  ! CVELZ will be multiplied by DTM in TESTCFL
138  current_state%local_cvel_z=max(current_state%local_cvel_z, &
139  current_state%abswmax(k)*current_state%global_grid%configuration%vertical%rdzn(k+1))
140  end do
141 #else
142  current_state%local_cvel_z=0.0_default_precision
143 #endif
144  call get_global_values(current_state%local_zumin, current_state%local_zumax, current_state%local_zvmin, &
145  current_state%local_zvmax, current_state%local_cvel_z, current_state%cvis, &
146  global_zumin, global_zumax, global_zvmin, global_zvmax, global_cvel_z, global_cvis, current_state%parallel)
147 
148  if (current_state%galilean_transformation) then
149  if (.not.current_state%fix_ugal)current_state%ugal=0.5_default_precision*(global_zumin+global_zumax)
150  if (.not.current_state%fix_vgal)current_state%vgal=0.5_default_precision*(global_zvmin+global_zvmax)
151  else
152  current_state%ugal=0.0_default_precision
153  current_state%vgal=0.0_default_precision
154  end if
155  current_state%cvel_z=global_cvel_z
156  current_state%cvel_x=max(abs(global_zumax-current_state%ugal), abs(global_zumin-current_state%ugal))
157  current_state%cvel_y=max(abs(global_zvmax-current_state%vgal), abs(global_zvmin-current_state%vgal))
158  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 55 of file cfltest.F90.

55  type(model_state_type), intent(inout), target :: current_state
56 
57  real(kind=DEFAULT_PRECISION) :: cfl_number
58 
59  if (mod(current_state%timestep, current_state%cfl_frequency) == 1 .or. &
60  current_state%timestep-current_state%start_timestep .le. current_state%cfl_frequency) then
61  current_state%cvel=0.0_default_precision
62  current_state%cvel_x=0.0_default_precision
63  current_state%cvel_y=0.0_default_precision
64  current_state%cvel_z=0.0_default_precision
65 
66  call perform_cfl_and_galilean_transformation_calculation(current_state)
67 
68  current_state%cvel=(current_state%cvel_x*current_state%global_grid%configuration%horizontal%cx+current_state%cvel_y*&
69  current_state%global_grid%configuration%horizontal%cy+current_state%cvel_z)*current_state%dtm
70  current_state%cvis=current_state%cvis*(current_state%dtm * 4)
71 
72  cfl_number=current_state%cvis/cvismax+current_state%cvel/cvelmax
73 
74  current_state%absolute_new_dtm=current_state%dtm
75  current_state%update_dtm=.false.
76  if (cfl_number .gt. 0.0_default_precision) then
77  if (cfl_number .lt. (1.0_default_precision-tollerance) .or. cfl_number .gt. (1.0_default_precision+tollerance)) then
78  current_state%absolute_new_dtm=current_state%dtm/cfl_number
79  end if
80  end if
81  end if
82  call update_dtm_based_on_absolute(current_state)
83  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)
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 91 of file cfltest.F90.

91  type(model_state_type), intent(inout), target :: current_state
92 
93  if (current_state%dtm .ne. current_state%absolute_new_dtm .and. &
94  (current_state%dtm .ne. dtmmax .or. current_state%absolute_new_dtm .lt. dtmmax)) then
95  current_state%update_dtm=.true.
96  current_state%dtm_new=min(current_state%dtm*(1.0_default_precision+rincmax), current_state%absolute_new_dtm, dtmmax)
97  if (current_state%parallel%my_rank==0) then
98  if (log_get_logging_level() .eq. log_debug) then
99  call log_log(log_debug, "dtm changed from "//trim(conv_to_string(current_state%dtm, 5))//" to "//&
100  trim(conv_to_string(current_state%dtm_new, 5)))
101  end if
102  if (current_state%dtm_new .lt. dtmmin) then
103  call log_log(log_error, "Timestep too small, dtmnew="//trim(conv_to_string(current_state%dtm_new, 5))//&
104  " dtmmin="//trim(conv_to_string(dtmmin, 5)))
105  end if
106  end if
107  end if
Here is the caller graph for this function:

Variable Documentation

◆ cvelmax

real(kind=default_precision) cfltest_mod::cvelmax
private

Definition at line 22 of file cfltest.F90.

◆ cvismax

real(kind=default_precision) cfltest_mod::cvismax
private

Definition at line 22 of file cfltest.F90.

◆ dtmmax

real(kind=default_precision) cfltest_mod::dtmmax
private

Definition at line 22 of file cfltest.F90.

◆ dtmmin

real(kind=default_precision) cfltest_mod::dtmmin
private

Definition at line 22 of file cfltest.F90.

◆ rincmax

real(kind=default_precision) cfltest_mod::rincmax
private

Definition at line 22 of file cfltest.F90.

◆ tollerance

real(kind=default_precision) cfltest_mod::tollerance
private

Definition at line 22 of file cfltest.F90.

22  real(kind=DEFAULT_PRECISION) :: tollerance, cvismax, cvelmax, dtmmax, dtmmin, rincmax