MONC
Functions/Subroutines | Variables
coriolis_mod Module Reference

This calculates the coriolis and mean pressure gradient terms which impact su and sv fields. More...

Functions/Subroutines

type(component_descriptor_type) function, public coriolis_get_descriptor ()
 Provides the descriptor back to the caller and is used in component registration. More...
 
subroutine initialisation_callback (current_state)
 Initialisation call back which will read in the coriolis configuration and set up the geostrophic winds. More...
 
subroutine timestep_callback (current_state)
 For each none halo cell this will calculate the coriolis terms for su and sv fields. More...
 

Variables

logical baroclinicity_use_geostrophic_shear
 
real(kind=default_precision), dimension(:), allocatable geostrophic_wind_x
 
real(kind=default_precision), dimension(:), allocatable geostrophic_wind_y
 
real(kind=default_precision) fcoriol
 

Detailed Description

This calculates the coriolis and mean pressure gradient terms which impact su and sv fields.

Function/Subroutine Documentation

◆ coriolis_get_descriptor()

type(component_descriptor_type) function, public coriolis_mod::coriolis_get_descriptor ( )

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

Returns
The termination check component descriptor

Definition at line 25 of file coriolis.F90.

25  coriolis_get_descriptor%name="coriolis"
26  coriolis_get_descriptor%version=0.1
27  coriolis_get_descriptor%initialisation=>initialisation_callback
28  coriolis_get_descriptor%timestep=>timestep_callback
Here is the call graph for this function:

◆ initialisation_callback()

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

Initialisation call back which will read in the coriolis configuration and set up the geostrophic winds.

Parameters
current_stateThe current model state

Definition at line 34 of file coriolis.F90.

34  type(model_state_type), target, intent(inout) :: current_state
35 
36  integer :: k
37 
38  baroclinicity_use_geostrophic_shear=options_get_logical(current_state%options_database, "baroclinicity_use_geostrophic_shear")
39  fcoriol=options_get_real(current_state%options_database, "fcoriol")
40  current_state%geostrophic_wind_rate_of_change_in_x=options_get_real(current_state%options_database, &
41  "geostrophic_wind_rate_of_change_in_x")
42  current_state%geostrophic_wind_rate_of_change_in_y=options_get_real(current_state%options_database, &
43  "geostrophic_wind_rate_of_change_in_y")
44  current_state%surface_geostrophic_wind_x=options_get_real(current_state%options_database, "surface_geostrophic_wind_x")
45  current_state%surface_geostrophic_wind_y=options_get_real(current_state%options_database, "surface_geostrophic_wind_y")
46 
47  allocate(geostrophic_wind_x(current_state%local_grid%size(z_index)), &
48  geostrophic_wind_y(current_state%local_grid%size(z_index)))
49 
50  do k=1,current_state%local_grid%size(z_index)
51  geostrophic_wind_x(k)=current_state%surface_geostrophic_wind_x
52  geostrophic_wind_y(k)=current_state%surface_geostrophic_wind_y
53  if (baroclinicity_use_geostrophic_shear) then
54  geostrophic_wind_x(k)=geostrophic_wind_x(k)+current_state%geostrophic_wind_rate_of_change_in_x*&
55  current_state%global_grid%configuration%vertical%zn(k)
56  geostrophic_wind_y(k)=geostrophic_wind_y(k)+current_state%geostrophic_wind_rate_of_change_in_y*&
57  current_state%global_grid%configuration%vertical%zn(k)
58  end if
59  end do
Here is the caller graph for this function:

◆ timestep_callback()

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

For each none halo cell this will calculate the coriolis terms for su and sv fields.

Parameters
current_stateThe current model state

Definition at line 65 of file coriolis.F90.

65  type(model_state_type), target, intent(inout) :: current_state
66 
67  integer :: k
68 
69  if (current_state%halo_column) then
70  if (.not. ((current_state%column_local_y == current_state%local_grid%halo_size(y_index) .and. &
71  current_state%column_local_x .le. current_state%local_grid%local_domain_end_index(x_index) .and. &
72  current_state%column_local_x .ge. current_state%local_grid%local_domain_start_index(x_index)-1) .or. &
73  (current_state%column_local_x == current_state%local_grid%halo_size(x_index) .and. &
74  current_state%column_local_y .ge. current_state%local_grid%local_domain_start_index(y_index) &
75  .and. current_state%column_local_y .le. current_state%local_grid%local_domain_end_index(y_index)) )) return
76  end if
77 
78  do k=2,current_state%local_grid%size(z_index)
79 #if defined(U_ACTIVE) && defined(V_ACTIVE)
80  current_state%su%data(k, current_state%column_local_y, current_state%column_local_x)=&
81  current_state%su%data(k, current_state%column_local_y, current_state%column_local_x)+fcoriol*&
82  (0.25_default_precision*(current_state%v%data(k, current_state%column_local_y, current_state%column_local_x)+&
83  current_state%v%data(k, current_state%column_local_y, current_state%column_local_x+1)+&
84  current_state%v%data(k, current_state%column_local_y-1, current_state%column_local_x)+&
85  current_state%v%data(k, current_state%column_local_y-1, current_state%column_local_x+1))+current_state%vgal-&
86  geostrophic_wind_y(k))
87 
88  current_state%sv%data(k, current_state%column_local_y, current_state%column_local_x)=&
89  current_state%sv%data(k, current_state%column_local_y, current_state%column_local_x)-fcoriol*&
90  (0.25_default_precision*(current_state%u%data(k, current_state%column_local_y, current_state%column_local_x)+&
91  current_state%u%data(k, current_state%column_local_y, current_state%column_local_x-1)+&
92  current_state%u%data(k, current_state%column_local_y+1, current_state%column_local_x)+&
93  current_state%u%data(k, current_state%column_local_y+1, current_state%column_local_x-1))+current_state%ugal-&
94  geostrophic_wind_x(k))
95 #endif
96  end do
Here is the caller graph for this function:

Variable Documentation

◆ baroclinicity_use_geostrophic_shear

logical coriolis_mod::baroclinicity_use_geostrophic_shear
private

Definition at line 14 of file coriolis.F90.

14  logical :: baroclinicity_use_geostrophic_shear

◆ fcoriol

real(kind=default_precision) coriolis_mod::fcoriol
private

Definition at line 16 of file coriolis.F90.

16  real(kind=DEFAULT_PRECISION) :: fcoriol

◆ geostrophic_wind_x

real(kind=default_precision), dimension(:), allocatable coriolis_mod::geostrophic_wind_x
private

Definition at line 15 of file coriolis.F90.

15  real(kind=DEFAULT_PRECISION), dimension(:), allocatable :: geostrophic_wind_x, geostrophic_wind_y

◆ geostrophic_wind_y

real(kind=default_precision), dimension(:), allocatable coriolis_mod::geostrophic_wind_y
private

Definition at line 15 of file coriolis.F90.