MONC
set_consistent_lowbc.F90
Go to the documentation of this file.
1 
10  use grids_mod, only : z_index
12 
14 
15  integer :: iqv ! index for vapour
16 
17  public set_consistent_lowbc_descriptor
18 
19 contains
23  set_consistent_lowbc_get_descriptor%name="set_consistent_lowbc"
28 
29  subroutine initialisation_callback(current_state)
30  ! copy of the initialisation callback from pwadvection, used to
31  ! identify whether acting on flow, theta or q fields
32  type(model_state_type), target, intent(inout) :: current_state
33 
34  ! Determine vapour index
35  if (.not. current_state%passive_q) then
36  iqv = get_q_index(standard_q_names%VAPOUR, 'lowerbc')
37  endif
38 
39  end subroutine initialisation_callback
40 
41  subroutine timestep_callback(current_state)
42  type(model_state_type), target, intent(inout) :: current_state
43 
44  integer :: current_x_index, current_y_index
45 
46  current_x_index=current_state%column_local_x
47  current_y_index=current_state%column_local_y
48 
49  if (current_state%halo_column) return
50 
51  call set_flow_lowbc(current_state, current_x_index, current_y_index)
52  if (current_state%th%active) then
53  call set_th_lowbc(current_state, current_x_index, current_y_index)
54  endif
55  if (current_state%number_q_fields .gt. 0) then
56  call set_q_lowbc(current_state, current_x_index, current_y_index)
57  endif
58  end subroutine timestep_callback
59 
60  subroutine set_flow_lowbc(current_state, current_x_index, current_y_index)
61  type(model_state_type), target, intent(inout) :: current_state
62  integer, intent(in) :: current_x_index, current_y_index
63 
64  if (current_state%use_viscosity_and_diffusion .and. &
65  current_state%use_surface_boundary_conditions) then
66 #ifdef U_ACTIVE
67  current_state%su%data(1, current_y_index, current_x_index)= &
68  -current_state%su%data(2, current_y_index, current_x_index)
69 #endif
70 #ifdef V_ACTIVE
71  current_state%sv%data(1, current_y_index, current_x_index)= &
72  -current_state%sv%data(2, current_y_index, current_x_index)
73 #endif
74  else
75 #ifdef U_ACTIVE
76  current_state%su%data(1, current_y_index, current_x_index)= &
77  current_state%su%data(2, current_y_index, current_x_index)
78 #endif
79 #ifdef V_ACTIVE
80  current_state%sv%data(1, current_y_index, current_x_index)= &
81  current_state%sv%data(2, current_y_index, current_x_index)
82 #endif
83  endif
84 
85  end subroutine set_flow_lowbc
86 
87  subroutine set_th_lowbc(current_state, current_x_index, current_y_index)
88  type(model_state_type), target, intent(inout) :: current_state
89  integer, intent(in) :: current_x_index, current_y_index
90 
91  if (current_state%use_surface_boundary_conditions) then
92  if (current_state%type_of_surface_boundary_conditions == prescribed_surface_fluxes) then
93  current_state%sth%data(1, current_y_index, current_x_index)= &
94  current_state%sth%data(2, current_y_index, current_x_index)
95  else
96  current_state%sth%data(1, current_y_index, current_x_index)= &
97  -current_state%sth%data(2, current_y_index, current_x_index)
98  endif
99  endif
100 
101  end subroutine set_th_lowbc
102 
103  subroutine set_q_lowbc(current_state, current_x_index, current_y_index)
104  type(model_state_type), target, intent(inout) :: current_state
105  integer, intent(in) :: current_x_index, current_y_index
106 
107  integer :: i
108 
109  if (current_state%number_q_fields .gt. 0) then
110  do n=1,current_state%number_q_fields
111  current_state%sq(n)%data(1, current_y_index, current_x_index)= &
112  current_state%sq(n)%data(2,current_y_index, current_x_index)
113  enddo
114  endif
115  if (current_state%use_surface_boundary_conditions .and. &
116  current_state%type_of_surface_boundary_conditions == prescribed_surface_values) then
117  current_state%sq(iqv)%data(1, current_y_index, current_x_index)= &
118  -(current_state%sq(iqv)%data(2,current_y_index, current_x_index))
119  endif
120 
121  end subroutine set_q_lowbc
122 
123 end module set_consistent_lowbc_mod
type(component_descriptor_type) function set_consistent_lowbc_get_descriptor()
Provides the descriptor back to the caller and is used in component registration. ...
type(standard_q_names_type), public standard_q_names
Definition: q_indices.F90:59
subroutine timestep_callback(current_state)
integer, parameter, public default_precision
MPI communication type which we use for the prognostic and calculation data.
Definition: datadefn.F90:17
integer, parameter, public z_index
Grid index parameters.
Definition: grids.F90:14
subroutine initialisation_callback(current_state)
Contains common definitions for the data and datatypes used by MONC.
Definition: datadefn.F90:2
The ModelState which represents the current state of a run.
Definition: state.F90:39
integer, parameter, public prescribed_surface_fluxes
Definition: state.F90:15
This manages the Q variables and specifically the mapping between names and the index that they are s...
Definition: q_indices.F90:2
Interfaces and types that MONC components must specify.
subroutine set_q_lowbc(current_state, current_x_index, current_y_index)
integer, parameter, public prescribed_surface_values
Definition: state.F90:15
real(kind=default_precision) function, public options_get_real(options_database, key, index)
Retrieves a real value from the database that matches the provided key.
Functionality to support the different types of grid and abstraction between global grids and local o...
Definition: grids.F90:5
This component sets the source term for the lowest level (Level 1) so that, depending on surface cons...
Manages the options database. Contains administration functions and deduce runtime options from the c...
subroutine set_flow_lowbc(current_state, current_x_index, current_y_index)
subroutine set_th_lowbc(current_state, current_x_index, current_y_index)
The model state which represents the current state of a run.
Definition: state.F90:2
integer function, public get_q_index(name, assigning_component)
Add in a new entry into the register if the name does not already exist or return the index of the pr...
Definition: q_indices.F90:112