MONC
Functions/Subroutines | Variables
thadvection_mod Module Reference

Specific theta advection, which involves the vertical advection of reference state and advection of mean baroclinicity. More...

Functions/Subroutines

type(component_descriptor_type) function, public thadvection_get_descriptor ()
 Provides a description of this component for the core to register. More...
 
subroutine initialisation_callback (current_state)
 Initialisation callback to set up the variables and data needed by the component. More...
 
subroutine timestep_callback (current_state)
 Timestep callback, will call the two separate procedures to do their advection if needed. More...
 
subroutine vertical_advection_of_reference_state (current_state, local_y, local_x)
 Vertical advection of the reference state. It doesn't seem consistent to do the advection in this way if TVD advection of the deviation from the reference state has been selected. Separate vertical advection of the reference state was introduced to improve energy conservation when carrying out idealized gravity wave simulations in a deep, dry isothermal layer, for which the difference in potential temp between top and bottom was of order 100K. In less extreme cases the benefits are unlikely to be significant and with TVD advection energy conservation has been compromised so the best way forward might be to recombine the reference state into l_th. More...
 
subroutine advection_of_mean_baroclinicity (current_state, local_y, local_x)
 Performs advection of the mean baroclinicity if appropriate. More...
 

Variables

logical baroclinicity_use_geostrophic_shear
 
real(kind=default_precision) fcoriol
 
real(kind=default_precision) fcoriol_over_g
 
real(kind=default_precision) rate_change_geostrophic_wind_x
 
real(kind=default_precision) rate_change_geostrophic_wind_y
 
real(kind=default_precision) multiplicative_factor_x
 
real(kind=default_precision) multiplicative_factor_y
 

Detailed Description

Specific theta advection, which involves the vertical advection of reference state and advection of mean baroclinicity.

Function/Subroutine Documentation

◆ advection_of_mean_baroclinicity()

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

Performs advection of the mean baroclinicity if appropriate.

Parameters
current_stateThe current model state
local_yThe local Y of the column
local_xThe local X of the column

Definition at line 89 of file thadvection.F90.

89  type(model_state_type), target, intent(inout) :: current_state
90  integer, intent(in) :: local_y, local_x
91 
92  integer :: k
93 
94  if (baroclinicity_use_geostrophic_shear) then
95  if (current_state%passive_q) then
96  if (current_state%use_anelastic_equations) then
97  do k=2, current_state%local_grid%size(z_index)
98  current_state%sth%data(k, local_y, local_x)=current_state%sth%data(k, local_y, local_x)+&
99  current_state%global_grid%configuration%vertical%thref(k)*fcoriol_over_g*&
100  ((current_state%v%data(k, local_y, local_x) + current_state%vgal) * rate_change_geostrophic_wind_x-&
101  (current_state%u%data(k, local_y, local_x) + current_state%ugal) * rate_change_geostrophic_wind_y)
102  end do
103  else
104  do k=2, current_state%local_grid%size(z_index)
105  current_state%sth%data(k, local_y, local_x)=current_state%sth%data(k, local_y, local_x)+&
106  ((current_state%v%data(k, local_y, local_x) + current_state%vgal) * multiplicative_factor_x-&
107  (current_state%u%data(k, local_y, local_x) + current_state%ugal) * multiplicative_factor_y)
108  end do
109  end if
110  else
111  call log_master_log(log_error, "The combination if baroclinicity and active q is not yet allowed")
112  end if
113  end if
Here is the caller graph for this function:

◆ initialisation_callback()

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

Initialisation callback to set up the variables and data needed by the component.

Parameters
current_stateThe current model state

Definition at line 35 of file thadvection.F90.

35  type(model_state_type), target, intent(inout) :: current_state
36 
37  baroclinicity_use_geostrophic_shear=options_get_logical(current_state%options_database, "baroclinicity_use_geostrophic_shear")
38  fcoriol=options_get_real(current_state%options_database, "fcoriol")
39  rate_change_geostrophic_wind_x=options_get_real(current_state%options_database, "rate_change_geostrophic_wind_x")
40  rate_change_geostrophic_wind_y=options_get_real(current_state%options_database, "rate_change_geostrophic_wind_y")
41  fcoriol_over_g = fcoriol/g
42  multiplicative_factor_x=rate_change_geostrophic_wind_x*current_state%thref0*fcoriol_over_g
43  multiplicative_factor_y=rate_change_geostrophic_wind_y*current_state%thref0*fcoriol_over_g
Here is the caller graph for this function:

◆ thadvection_get_descriptor()

type(component_descriptor_type) function, public thadvection_mod::thadvection_get_descriptor ( )

Provides a description of this component for the core to register.

Returns
The descriptor containing registration information for this component

Definition at line 26 of file thadvection.F90.

26  thadvection_get_descriptor%name="th_advection"
27  thadvection_get_descriptor%version=0.1
28  thadvection_get_descriptor%initialisation=>initialisation_callback
29  thadvection_get_descriptor%timestep=>timestep_callback
Here is the call graph for this function:

◆ timestep_callback()

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

Timestep callback, will call the two separate procedures to do their advection if needed.

Parameters
current_stateThe current model state

Definition at line 49 of file thadvection.F90.

49  type(model_state_type), target, intent(inout) :: current_state
50 
51  call vertical_advection_of_reference_state(current_state, current_state%column_local_y, current_state%column_local_x)
52  call advection_of_mean_baroclinicity(current_state, current_state%column_local_y, current_state%column_local_x)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ vertical_advection_of_reference_state()

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

Vertical advection of the reference state. It doesn't seem consistent to do the advection in this way if TVD advection of the deviation from the reference state has been selected. Separate vertical advection of the reference state was introduced to improve energy conservation when carrying out idealized gravity wave simulations in a deep, dry isothermal layer, for which the difference in potential temp between top and bottom was of order 100K. In less extreme cases the benefits are unlikely to be significant and with TVD advection energy conservation has been compromised so the best way forward might be to recombine the reference state into l_th.

Parameters
current_stateThe current model state
local_yThe local Y of the column
local_xThe local X of the column

Definition at line 65 of file thadvection.F90.

65  type(model_state_type), target, intent(inout) :: current_state
66  integer, intent(in) :: local_y, local_x
67 
68  integer :: k
69  real(kind=DEFAULT_PRECISION) :: sctmp1, sctmp2
70 
71  if (current_state%use_anelastic_equations) then
72  ! This code only needs to be executed if anelastic, otherwise THREF is constant and the increment calculated here is zero
73  do k=2, current_state%local_grid%size(z_index)
74  sctmp1=current_state%global_grid%configuration%vertical%tzc1(k)*2.0_default_precision*&
75  current_state%global_grid%configuration%vertical%dthref(k-1)
76  sctmp2=current_state%global_grid%configuration%vertical%tzc2(k)*2.0_default_precision*&
77  current_state%global_grid%configuration%vertical%dthref(k)
78  current_state%sth%data(k, local_y, local_x)=current_state%sth%data(k, local_y, local_x)-(sctmp1*&
79  current_state%w%data(k-1, local_y, local_x) + sctmp2*current_state%w%data(k, local_y, local_x))
80  end do
81  end if
Here is the caller graph for this function:

Variable Documentation

◆ baroclinicity_use_geostrophic_shear

logical thadvection_mod::baroclinicity_use_geostrophic_shear
private

Definition at line 16 of file thadvection.F90.

16  logical :: baroclinicity_use_geostrophic_shear

◆ fcoriol

real(kind=default_precision) thadvection_mod::fcoriol
private

Definition at line 17 of file thadvection.F90.

17  real(kind=DEFAULT_PRECISION) :: fcoriol, fcoriol_over_g, rate_change_geostrophic_wind_x, rate_change_geostrophic_wind_y, &
18  multiplicative_factor_x, multiplicative_factor_y

◆ fcoriol_over_g

real(kind=default_precision) thadvection_mod::fcoriol_over_g
private

Definition at line 17 of file thadvection.F90.

◆ multiplicative_factor_x

real(kind=default_precision) thadvection_mod::multiplicative_factor_x
private

Definition at line 17 of file thadvection.F90.

◆ multiplicative_factor_y

real(kind=default_precision) thadvection_mod::multiplicative_factor_y
private

Definition at line 17 of file thadvection.F90.

◆ rate_change_geostrophic_wind_x

real(kind=default_precision) thadvection_mod::rate_change_geostrophic_wind_x
private

Definition at line 17 of file thadvection.F90.

◆ rate_change_geostrophic_wind_y

real(kind=default_precision) thadvection_mod::rate_change_geostrophic_wind_y
private

Definition at line 17 of file thadvection.F90.