MONC
Functions/Subroutines | Variables
buoyancy_mod Module Reference

Calculates buoyancy terms for the SW field. More...

Functions/Subroutines

type(component_descriptor_type) function, public buoyancy_get_descriptor ()
 Provides the descriptor back to the caller and is used in component registration. More...
 
subroutine field_information_retrieval_callback (current_state, name, field_information)
 Field information retrieval callback, this returns information for a specific components published field. More...
 
subroutine field_value_retrieval_callback (current_state, name, field_value)
 Field value retrieval callback, this returns the value of a specific published field. More...
 
subroutine initialisation_callback (current_state)
 The initialisation callback sets up the buoyancy coefficient. More...
 
subroutine finalisation_callback (current_state)
 
subroutine timestep_callback (current_state)
 Called for each column per timestep this will calculate the buoyancy terms for the SW field. More...
 
subroutine save_precomponent_tendencies (current_state, cxn, cyn, txn, tyn)
 Save the 3d tendencies coming into the component. More...
 
subroutine compute_component_tendencies (current_state, cxn, cyn, txn, tyn)
 Computation of component tendencies. More...
 
subroutine set_published_field_value (field_value, real_1d_field, real_2d_field, real_3d_field)
 Sets the published field value from the temporary diagnostic values held by this component. More...
 

Variables

real(kind=default_precision), dimension(:), allocatable w_buoyancy
 
real(kind=default_precision) g_over_2
 
integer iqv
 
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_w
 
logical l_tend_3d_w
 
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_w
 
logical l_tend_pr_tot_w
 
integer diagnostic_generation_frequency
 

Detailed Description

Calculates buoyancy terms for the SW field.

Function/Subroutine Documentation

◆ buoyancy_get_descriptor()

type(component_descriptor_type) function, public buoyancy_mod::buoyancy_get_descriptor

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

Returns
The termination check component descriptor

Definition at line 40 of file buoyancy.F90.

41  buoyancy_get_descriptor%name="buoyancy"
42  buoyancy_get_descriptor%version=0.1
43  buoyancy_get_descriptor%initialisation=>initialisation_callback
44  buoyancy_get_descriptor%timestep=>timestep_callback
45  buoyancy_get_descriptor%finalisation=>finalisation_callback
46 
47  buoyancy_get_descriptor%field_value_retrieval=>field_value_retrieval_callback
48  buoyancy_get_descriptor%field_information_retrieval=>field_information_retrieval_callback
49  allocate(buoyancy_get_descriptor%published_fields(1+1+1))
50 
51  buoyancy_get_descriptor%published_fields(1)="w_buoyancy"
52 
53  buoyancy_get_descriptor%published_fields(1+1)= "tend_w_buoyancy_3d_local"
54 
55  buoyancy_get_descriptor%published_fields(1+1+1)= "tend_w_buoyancy_profile_total_local"
56 
Here is the call graph for this function:

◆ compute_component_tendencies()

subroutine buoyancy_mod::compute_component_tendencies ( type(model_state_type), intent(inout), target  current_state,
integer, intent(in)  cxn,
integer, intent(in)  cyn,
integer, intent(in)  txn,
integer, intent(in)  tyn 
)
private

Computation of component tendencies.

Parameters
current_stateCurrent model state
cxnThe current slice, x, index
cynThe current column, y, index.
txntarget_x_index
tyntarget_y_index

Definition at line 286 of file buoyancy.F90.

287  type(model_state_type), target, intent(inout) :: current_state
288  integer, intent(in) :: cxn, cyn, txn, tyn
289 
290  ! Calculate change in tendency due to component
291  if (l_tend_3d_w) then
292  tend_3d_w(:,tyn,txn)=current_state%sw%data(:,cyn,cxn) - tend_3d_w(:,tyn,txn)
293  endif
294 
295  ! Add local tendency fields to the profile total
296  if (l_tend_pr_tot_w) then
297  tend_pr_tot_w(:)=tend_pr_tot_w(:) + tend_3d_w(:,tyn,txn)
298  endif
299 
Here is the caller graph for this function:

◆ field_information_retrieval_callback()

subroutine buoyancy_mod::field_information_retrieval_callback ( type(model_state_type), intent(inout), target  current_state,
character(len=*), intent(in)  name,
type(component_field_information_type), intent(out)  field_information 
)
private

Field information retrieval callback, this returns information for a specific components published field.

Parameters
current_stateCurrent model state
nameThe name of the field to retrieve information for
field_informationPopulated with information about the field

Definition at line 63 of file buoyancy.F90.

64  type(model_state_type), target, intent(inout) :: current_state
65  character(len=*), intent(in) :: name
66  type(component_field_information_type), intent(out) :: field_information
67  integer :: strcomp
68 
69  ! Field description is the same regardless of the specific field being retrieved
70  field_information%field_type=component_array_field_type
71  field_information%data_type=component_double_data_type
72  field_information%number_dimensions=1
73  field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
74  field_information%enabled=.true.
75 
76  ! Field information for 3d
77  strcomp=index(name, "buoyancy_3d_local")
78  if (strcomp .ne. 0) then
79  field_information%field_type=component_array_field_type
80  field_information%number_dimensions=3
81  field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
82  field_information%dimension_sizes(2)=current_state%local_grid%size(y_index)
83  field_information%dimension_sizes(3)=current_state%local_grid%size(x_index)
84  field_information%data_type=component_double_data_type
85 
86  if (name .eq. "tend_w_buoyancy_3d_local") then
87  field_information%enabled=l_tend_3d_w
88  else
89  field_information%enabled=.true.
90  end if
91 
92  end if !end 3d check
93 
94  ! Field information for profiles
95  strcomp=index(name, "buoyancy_profile_total_local")
96  if (strcomp .ne. 0) then
97  field_information%field_type=component_array_field_type
98  field_information%number_dimensions=1
99  field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
100  field_information%data_type=component_double_data_type
101 
102  if (name .eq. "tend_w_buoyancy_profile_total_local") then
103  field_information%enabled=l_tend_pr_tot_w
104  else
105  field_information%enabled=.true.
106  end if
107 
108  end if !end profile check
109 
Here is the caller graph for this function:

◆ field_value_retrieval_callback()

subroutine buoyancy_mod::field_value_retrieval_callback ( type(model_state_type), intent(inout), target  current_state,
character(len=*), intent(in)  name,
type(component_field_value_type), intent(out)  field_value 
)
private

Field value retrieval callback, this returns the value of a specific published field.

Parameters
current_stateCurrent model state
nameThe name of the field to retrieve the value for
field_valuePopulated with the value of the field

Definition at line 117 of file buoyancy.F90.

118  type(model_state_type), target, intent(inout) :: current_state
119  character(len=*), intent(in) :: name
120  type(component_field_value_type), intent(out) :: field_value
121 
122  if (name .eq. "w_buoyancy") then
123  allocate(field_value%real_1d_array(size(w_buoyancy)), source=w_buoyancy)
124 
125  ! 3d Tendency Fields
126  else if (name .eq. "tend_w_buoyancy_3d_local" .and. allocated(tend_3d_w)) then
127  call set_published_field_value(field_value, real_3d_field=tend_3d_w)
128 
129  ! Profile Tendency Fields
130  else if (name .eq. "tend_w_buoyancy_profile_total_local" .and. allocated(tend_pr_tot_w)) then
131  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_w)
132  end if
133 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ finalisation_callback()

subroutine buoyancy_mod::finalisation_callback ( type(model_state_type), intent(inout), target  current_state)
private

Definition at line 179 of file buoyancy.F90.

180  type(model_state_type), target, intent(inout) :: current_state
181 
182  if (allocated(w_buoyancy)) deallocate(w_buoyancy)
183 
184  if (allocated(tend_3d_w)) deallocate(tend_3d_w)
185 
186  if (allocated(tend_pr_tot_w)) deallocate(tend_pr_tot_w)
187 
Here is the caller graph for this function:

◆ initialisation_callback()

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

The initialisation callback sets up the buoyancy coefficient.

Parameters
current_stateThe current model state

Definition at line 139 of file buoyancy.F90.

140  type(model_state_type), target, intent(inout) :: current_state
141 
142  integer :: z_size
143 
144  if (.not. current_state%passive_q .and. current_state%number_q_fields > 0)then
145  if (.not. allocated(current_state%cq))then
146  allocate(current_state%cq(current_state%number_q_fields))
147  current_state%cq=0.0_default_precision
148  end if
149  iqv = get_q_index(standard_q_names%VAPOUR, 'buoyancy')
150  current_state%cq(iqv) = ratio_mol_wts-1.0
151  end if
152 
153  g_over_2 = 0.5_default_precision*g
154  z_size=current_state%global_grid%size(z_index)
155  allocate(w_buoyancy(z_size))
156 
157  ! Tendency logicals
158  l_tend_pr_tot_w = current_state%w%active
159  l_tend_3d_w = current_state%w%active .or. l_tend_pr_tot_w
160 
161  ! Allocate 3d tendency fields upon availability
162  if (l_tend_3d_w) then
163  allocate( tend_3d_w(current_state%local_grid%size(z_index), &
164  current_state%local_grid%size(y_index), &
165  current_state%local_grid%size(x_index) ) )
166  endif
167 
168  ! Allocate profile tendency fields upon availability
169  if (l_tend_pr_tot_w) then
170  allocate( tend_pr_tot_w(current_state%local_grid%size(z_index)) )
171  endif
172 
173  ! Save the sampling_frequency to force diagnostic calculation on select time steps
174  diagnostic_generation_frequency=options_get_integer(current_state%options_database, "sampling_frequency")
175 
Here is the caller graph for this function:

◆ save_precomponent_tendencies()

subroutine buoyancy_mod::save_precomponent_tendencies ( type(model_state_type), intent(in), target  current_state,
integer, intent(in)  cxn,
integer, intent(in)  cyn,
integer, intent(in)  txn,
integer, intent(in)  tyn 
)
private

Save the 3d tendencies coming into the component.

Parameters
current_stateCurrent model state
cxnThe current slice, x, index
cynThe current column, y, index.
txntarget_x_index
tyntarget_y_index

Definition at line 268 of file buoyancy.F90.

269  type(model_state_type), target, intent(in) :: current_state
270  integer, intent(in) :: cxn, cyn, txn, tyn
271 
272  ! Save 3d tendency fields upon request (of 3d or profiles) and availability
273  if (l_tend_3d_w) then
274  tend_3d_w(:,tyn,txn)=current_state%sw%data(:,cyn,cxn)
275  endif
276 
Here is the caller graph for this function:

◆ set_published_field_value()

subroutine buoyancy_mod::set_published_field_value ( type(component_field_value_type), intent(inout)  field_value,
real(kind=default_precision), dimension(:), optional  real_1d_field,
real(kind=default_precision), dimension(:,:), optional  real_2d_field,
real(kind=default_precision), dimension(:,:,:), optional  real_3d_field 
)
private

Sets the published field value from the temporary diagnostic values held by this component.

Parameters
field_valuePopulated with the value of the field
real_1d_fieldOptional one dimensional real of values to publish
real_2d_fieldOptional two dimensional real of values to publish

Definition at line 307 of file buoyancy.F90.

308  type(component_field_value_type), intent(inout) :: field_value
309  real(kind=default_precision), dimension(:), optional :: real_1d_field
310  real(kind=default_precision), dimension(:,:), optional :: real_2d_field
311  real(kind=default_precision), dimension(:,:,:), optional :: real_3d_field
312 
313  if (present(real_1d_field)) then
314  allocate(field_value%real_1d_array(size(real_1d_field)), source=real_1d_field)
315  else if (present(real_2d_field)) then
316  allocate(field_value%real_2d_array(size(real_2d_field, 1), size(real_2d_field, 2)), source=real_2d_field)
317  else if (present(real_3d_field)) then
318  allocate(field_value%real_3d_array(size(real_3d_field, 1), size(real_3d_field, 2), size(real_3d_field, 3)), &
319  source=real_3d_field)
320  end if
Here is the caller graph for this function:

◆ timestep_callback()

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

Called for each column per timestep this will calculate the buoyancy terms for the SW field.

Parameters
current_stateThe current model state
target_(x/y)_indexThis is the index with the halos subtracted. This is needed so that diagnostic does not include halos and to prevent array out-of-bounds

Definition at line 195 of file buoyancy.F90.

196  type(model_state_type), target, intent(inout) :: current_state
197 
198  integer :: k, n
199  integer :: current_x_index, current_y_index, target_x_index, target_y_index
200 
201  current_x_index=current_state%column_local_x
202  current_y_index=current_state%column_local_y
203  target_y_index=current_y_index-current_state%local_grid%halo_size(y_index)
204  target_x_index=current_x_index-current_state%local_grid%halo_size(x_index)
205 
206  ! Zero profile tendency totals on first instance in the sum
207  if (current_state%first_timestep_column) then
208  if (l_tend_pr_tot_w) then
209  tend_pr_tot_w(:)=0.0_default_precision
210  endif
211  endif ! zero totals
212 
213  if (mod(current_state%timestep, diagnostic_generation_frequency) == 0 .and. .not. current_state%halo_column) then
214  call save_precomponent_tendencies(current_state, current_x_index, current_y_index, target_x_index, target_y_index)
215  end if
216 
217 
218 #ifdef W_ACTIVE
219  if (.not. current_state%passive_th .and. current_state%th%active) then
220  do k=2,current_state%local_grid%size(z_index)-1
221  w_buoyancy(k)=(0.5_default_precision*current_state%global_grid%configuration%vertical%buoy_co(k))*&
222  (current_state%th%data(k, current_state%column_local_y, current_state%column_local_x)&
223  +current_state%th%data(k+1, current_state%column_local_y, current_state%column_local_x))
224  current_state%sw%data(k, current_state%column_local_y, current_state%column_local_x)=&
225  current_state%sw%data(k, current_state%column_local_y, current_state%column_local_x)+w_buoyancy(k)
226  end do
227  end if
228  if (.not. current_state%passive_q .and. current_state%number_q_fields .gt. 0) then
229  if (current_state%use_anelastic_equations) then
230  do n=1,current_state%number_q_fields
231  do k=2,current_state%local_grid%size(z_index)-1
232  current_state%sw%data(k, current_state%column_local_y, current_state%column_local_x)=&
233  current_state%sw%data(k, current_state%column_local_y, current_state%column_local_x)+&
234  (0.5_default_precision*current_state%global_grid%configuration%vertical%buoy_co(k))*&
235  current_state%cq(n)* (current_state%global_grid%configuration%vertical%thref(k)*&
236  current_state%q(n)%data(k, current_state%column_local_y, current_state%column_local_x)+&
237  current_state%global_grid%configuration%vertical%thref(k+1)*&
238  current_state%q(n)%data(k+1, current_state%column_local_y, current_state%column_local_x))
239  end do
240  end do
241  else
242  do n=1,current_state%number_q_fields
243  do k=2,current_state%local_grid%size(z_index)-1
244  current_state%sw%data(k, current_state%column_local_y, current_state%column_local_x)=&
245  current_state%sw%data(k, current_state%column_local_y, current_state%column_local_x)+&
246  g_over_2*current_state%cq(n)*&
247  (current_state%q(n)%data(k, current_state%column_local_y, current_state%column_local_x)+&
248  current_state%q(n)%data(k+1, current_state%column_local_y, current_state%column_local_x))
249  end do
250  end do
251  end if
252  end if
253 #endif
254 
255  if (mod(current_state%timestep, diagnostic_generation_frequency) == 0 .and. .not. current_state%halo_column) then
256  call compute_component_tendencies(current_state, current_x_index, current_y_index, target_x_index, target_y_index)
257  end if
258 
Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ diagnostic_generation_frequency

integer buoyancy_mod::diagnostic_generation_frequency
private

Definition at line 32 of file buoyancy.F90.

32  integer :: diagnostic_generation_frequency

◆ g_over_2

real(kind=default_precision) buoyancy_mod::g_over_2
private

Definition at line 20 of file buoyancy.F90.

20  real(kind=default_precision) :: g_over_2

◆ iqv

integer buoyancy_mod::iqv
private

Definition at line 22 of file buoyancy.F90.

22  integer :: iqv ! Index for water vapour

◆ l_tend_3d_w

logical buoyancy_mod::l_tend_3d_w
private

Definition at line 27 of file buoyancy.F90.

27  logical :: l_tend_3d_w

◆ l_tend_pr_tot_w

logical buoyancy_mod::l_tend_pr_tot_w
private

Definition at line 30 of file buoyancy.F90.

30  logical :: l_tend_pr_tot_w

◆ tend_3d_w

real(kind=default_precision), dimension(:,:,:), allocatable buoyancy_mod::tend_3d_w
private

Definition at line 26 of file buoyancy.F90.

26  real(kind=default_precision), dimension(:,:,:), allocatable :: tend_3d_w

◆ tend_pr_tot_w

real(kind=default_precision), dimension(:), allocatable buoyancy_mod::tend_pr_tot_w
private

Definition at line 29 of file buoyancy.F90.

29  real(kind=default_precision), dimension(:), allocatable :: tend_pr_tot_w

◆ w_buoyancy

real(kind=default_precision), dimension(:), allocatable buoyancy_mod::w_buoyancy
private

Definition at line 18 of file buoyancy.F90.

18  real(kind=default_precision), dimension(:), allocatable :: w_buoyancy
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