MONC
Functions/Subroutines | Variables
lowerbc_mod Module Reference

This sets the lower boundary conditions for theta and the q variables. More...

Functions/Subroutines

type(component_descriptor_type) function, public lowerbc_get_descriptor ()
 Descriptor of this component for registration. More...
 
subroutine initialisation_callback (current_state)
 
subroutine finalisation_callback (current_state)
 
subroutine allocate_applicable_fields (current_state)
 
subroutine timestep_callback (current_state)
 
subroutine compute_lower_boundary_conditions (current_state, current_y_index, current_x_index, zu, zv, zth, th, zq, q)
 
subroutine register_async_wrapping_recv_requests (current_state)
 Registers asynchronous wrapping recv requests as needed. More...
 
subroutine complete_async_wrapping (current_state, zth, zq)
 Completes the asynchronous wrapping if required for periodic boundary conditions. More...
 
subroutine package_y_wrapping_send_buffer (current_state, zth, zq, first_y_index, second_y_index)
 Packages theta and Q fields (if enabled) into the send buffer for Y. More...
 
subroutine package_x_wrapping_send_buffer (current_state, zth, zq, first_x_index, second_x_index)
 Packages theta and Q fields (if enabled) into the send buffer for X. More...
 
subroutine unpackage_y_wrapping_recv_buffer (current_state, zth, zq, first_y_index, second_y_index)
 Unpackages theta and Q fields from the receive buffer into the fields themselves (if enabled) for Y. More...
 
subroutine unpackage_x_wrapping_recv_buffer (current_state, zth, zq, first_x_index, second_x_index)
 Unpackages theta and Q fields from the receive buffer into the fields themselves (if enabled) for X. More...
 
subroutine handle_convective_fluxes (current_state, current_y_index, current_x_index, horizontal_velocity_at_k2, th, q)
 
real(kind=default_precision) function look (current_state, vel)
 
subroutine handle_neutral_fluxes (current_state, current_y_index, current_x_index, horizontal_velocity_at_k2, th, q)
 
subroutine handle_stable_fluxes (current_state, current_y_index, current_x_index, horizontal_velocity_at_k2, th, q)
 
subroutine compute_using_fixed_surface_fluxes (current_state, current_y_index, current_x_index, horizontal_velocity_at_k2, th, q)
 
subroutine compute_using_fixed_surface_temperature (current_state, current_y_index, current_x_index, horizontal_velocity_at_k2, zth, th, zq, q)
 
subroutine simple_boundary_values (current_state, current_y_index, current_x_index, th, q)
 
integer function mostbc (current_state, delu, delt, z1, ustrdg, tstrdg)
 Solves the Monin-Obukhov equations in the case of specified surface values of temperature and mixing ratio,combined into a specified value of virtual temperature. It is a modified version of the subroutine described in Bull and Derbyshire (TDN197) based on the assumption that the similarity functions and roughness lengths for temperature and mixing ratio are the same. In that case, all the original theory can be used if we replace temperature by virtual temperature. The form of the non-dimensionalised wind shear used is 1.0 + BETAM z/L for the stable case and (1.0 - GAMMAM z/L)**1/4 for the unstable case. The form of the non-dimensionalised temperature gradient used is 1.0 + BETAH z/L for the stable case and (1.0 - GAMMAH z/L)**1/2 for the unstable case. More...
 
integer function solve_monin_obukhov_unstable_case (delu, delt, ellmocon, ustrdg, tstrdg, vertical_grid)
 
integer function solve_monin_obukhov_stable_case (delu, delt, zlogm, cmbc, ustrdg, tstrdg)
 

Variables

integer, parameter convergence_success =1
 
integer, parameter convergence_richardson_too_large =2
 
integer, parameter convergence_failure =3
 
real(kind=default_precision), parameter smth = 0.05_DEFAULT_PRECISION
 
real(kind=default_precision), parameter tolm =1.0E-4_DEFAULT_PRECISION
 
real(kind=default_precision), parameter tolt =1.0E-4_DEFAULT_PRECISION
 
real(kind=default_precision) tstrcona
 
real(kind=default_precision) rhmbc
 
real(kind=default_precision) ddbc
 
real(kind=default_precision) ddbc_x4
 
real(kind=default_precision) eecon
 
real(kind=default_precision) r2ddbc
 
real(kind=default_precision) rcmbc
 
real(kind=default_precision) tstrconb
 
real(kind=default_precision) x4con
 
real(kind=default_precision) xx0con
 
real(kind=default_precision) y2con
 
real(kind=default_precision) yy0con
 
real(kind=default_precision) viscous_courant_coefficient
 
real(kind=default_precision), dimension(:,:,:), allocatable x_wrapping_send_buffer
 
real(kind=default_precision), dimension(:,:,:), allocatable y_wrapping_send_buffer
 
real(kind=default_precision), dimension(:,:,:), allocatable x_wrapping_recv_buffer
 
real(kind=default_precision), dimension(:,:,:), allocatable y_wrapping_recv_buffer
 
integer iqv
 
integer, dimension(4) wrapping_comm_requests
 
integer y_wrapping_target_id
 
integer x_wrapping_target_id
 

Detailed Description

This sets the lower boundary conditions for theta and the q variables.

Function/Subroutine Documentation

◆ allocate_applicable_fields()

subroutine lowerbc_mod::allocate_applicable_fields ( type(model_state_type), intent(inout), target  current_state)
private

Definition at line 157 of file lowerbc.F90.

157  type(model_state_type), target, intent(inout) :: current_state
158 
159  integer :: z_size, y_size, x_size, i
160 
161  z_size=current_state%local_grid%size(z_index) + current_state%local_grid%halo_size(z_index) * 2
162  y_size=current_state%local_grid%size(y_index) + current_state%local_grid%halo_size(y_index) * 2
163  x_size=current_state%local_grid%size(x_index) + current_state%local_grid%halo_size(x_index) * 2
164 
165  allocate(current_state%dis%data(z_size, y_size, x_size), &
166  current_state%dis_th%data(z_size, y_size, x_size), current_state%disq(current_state%number_q_fields))
167 
168  do i=1,current_state%number_q_fields
169  allocate(current_state%disq(i)%data(z_size, y_size, x_size))
170  end do
Here is the caller graph for this function:

◆ complete_async_wrapping()

subroutine lowerbc_mod::complete_async_wrapping ( type(model_state_type), intent(inout), target  current_state,
type(prognostic_field_type), intent(inout)  zth,
type(prognostic_field_type), dimension(:), intent(inout)  zq 
)
private

Completes the asynchronous wrapping if required for periodic boundary conditions.

Parameters
current_stateThe current model state
zthTemperature field
zqQ fields

Definition at line 283 of file lowerbc.F90.

283  type(model_state_type), target, intent(inout) :: current_state
284  type(prognostic_field_type), intent(inout) :: zth
285  type(prognostic_field_type), dimension(:), intent(inout) :: zq
286 
287  integer :: ierr, n
288 
289  if (allocated(x_wrapping_send_buffer) .or. allocated(y_wrapping_send_buffer)) then
290  if (allocated(y_wrapping_send_buffer)) then
291  if (current_state%parallel%my_coords(y_index) == 0) then
292  call package_y_wrapping_send_buffer(current_state, zth, zq, current_state%local_grid%local_domain_start_index(y_index),&
293  current_state%local_grid%local_domain_start_index(y_index)+1)
294  else
295  call package_y_wrapping_send_buffer(current_state, zth, zq, current_state%local_grid%local_domain_end_index(y_index)-1,&
296  current_state%local_grid%local_domain_end_index(y_index))
297  end if
298  call mpi_isend(y_wrapping_send_buffer, size(y_wrapping_send_buffer), precision_type, &
299  y_wrapping_target_id, 0, current_state%parallel%neighbour_comm, &
300  wrapping_comm_requests(2), ierr)
301  end if
302  if (allocated(x_wrapping_send_buffer)) then
303  if (current_state%parallel%my_coords(x_index) == 0) then
304  call package_x_wrapping_send_buffer(current_state, zth, zq, current_state%local_grid%local_domain_start_index(x_index),&
305  current_state%local_grid%local_domain_start_index(x_index)+1)
306  else
307  call package_x_wrapping_send_buffer(current_state, zth, zq, current_state%local_grid%local_domain_end_index(x_index)-1,&
308  current_state%local_grid%local_domain_end_index(x_index))
309  end if
310  call mpi_isend(x_wrapping_send_buffer, size(x_wrapping_send_buffer), precision_type, &
311  x_wrapping_target_id, 0, current_state%parallel%neighbour_comm, &
312  wrapping_comm_requests(4), ierr)
313  end if
314 
315  ! If send buffer is allocated then recv buffer is allocated, therefore only test the send buffer here and assume recv
316  call mpi_waitall(4, wrapping_comm_requests, mpi_statuses_ignore, ierr)
317  wrapping_comm_requests=mpi_request_null
318  if (allocated(y_wrapping_recv_buffer)) then
319  if (current_state%parallel%my_coords(y_index) == 0) then
320  call unpackage_y_wrapping_recv_buffer(current_state, zth, zq, 1, 2)
321  else
322  call unpackage_y_wrapping_recv_buffer(current_state, zth, zq, &
323  current_state%local_grid%local_domain_end_index(y_index)+1, &
324  current_state%local_grid%local_domain_end_index(y_index)+2)
325  end if
326  end if
327  if (allocated(x_wrapping_recv_buffer)) then
328  if (current_state%parallel%my_coords(x_index) == 0) then
329  call unpackage_x_wrapping_recv_buffer(current_state, zth, zq, 1, 2)
330  else
331  call unpackage_x_wrapping_recv_buffer(current_state, zth, zq, &
332  current_state%local_grid%local_domain_end_index(x_index)+1, &
333  current_state%local_grid%local_domain_end_index(x_index)+2)
334  end if
335  end if
336  end if
337  if (current_state%parallel%my_rank == y_wrapping_target_id) then
338  if (current_state%th%active) then
339  zth%data(1,1,:)=zth%data(1, current_state%local_grid%local_domain_end_index(y_index)-1, :)
340  zth%data(1,2,:)=zth%data(1, current_state%local_grid%local_domain_end_index(y_index), :)
341  zth%data(1,current_state%local_grid%local_domain_end_index(y_index)+1,:)=&
342  zth%data(1, current_state%local_grid%local_domain_start_index(y_index), :)
343  zth%data(1,current_state%local_grid%local_domain_end_index(y_index)+2,:)=&
344  zth%data(1, current_state%local_grid%local_domain_start_index(y_index)+1, :)
345  end if
346  if (current_state%number_q_fields .gt. 0) then
347  do n=1, current_state%number_q_fields
348  zq(n)%data(1,1,:)=zq(n)%data(1, current_state%local_grid%local_domain_end_index(y_index)-1, :)
349  zq(n)%data(1,2,:)=zq(n)%data(1, current_state%local_grid%local_domain_end_index(y_index), :)
350  zq(n)%data(1,current_state%local_grid%local_domain_end_index(y_index)+1,:)=&
351  zq(n)%data(1, current_state%local_grid%local_domain_start_index(y_index), :)
352  zq(n)%data(1,current_state%local_grid%local_domain_end_index(y_index)+2,:)=&
353  zq(n)%data(1, current_state%local_grid%local_domain_start_index(y_index)+1, :)
354  end do
355  end if
356  end if
357 
358  if (current_state%parallel%my_rank == x_wrapping_target_id) then
359  if (current_state%th%active) then
360  zth%data(1,:,1)=zth%data(1,:,current_state%local_grid%local_domain_end_index(x_index)-1)
361  zth%data(1,:,2)=zth%data(1,:,current_state%local_grid%local_domain_end_index(x_index))
362  zth%data(1,:,current_state%local_grid%local_domain_end_index(x_index)+1)=&
363  zth%data(1,:,current_state%local_grid%local_domain_start_index(x_index))
364  zth%data(1,:,current_state%local_grid%local_domain_end_index(x_index)+2)=&
365  zth%data(1,:,current_state%local_grid%local_domain_start_index(x_index)+1)
366  end if
367  if (current_state%number_q_fields .gt. 0) then
368  do n=1, current_state%number_q_fields
369  zq(n)%data(1,:,1)=zq(n)%data(1,:,current_state%local_grid%local_domain_end_index(x_index)-1)
370  zq(n)%data(1,:,2)=zq(n)%data(1,:,current_state%local_grid%local_domain_end_index(x_index))
371  zq(n)%data(1,:,current_state%local_grid%local_domain_end_index(x_index)+1)=&
372  zq(n)%data(1,:,current_state%local_grid%local_domain_start_index(x_index))
373  zq(n)%data(1,:,current_state%local_grid%local_domain_end_index(x_index)+2)=&
374  zq(n)%data(1,:,current_state%local_grid%local_domain_start_index(x_index)+1)
375  end do
376  end if
377  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ compute_lower_boundary_conditions()

subroutine lowerbc_mod::compute_lower_boundary_conditions ( type(model_state_type), intent(inout), target  current_state,
integer, intent(in)  current_y_index,
integer, intent(in)  current_x_index,
type(prognostic_field_type), intent(inout)  zu,
type(prognostic_field_type), intent(inout)  zv,
type(prognostic_field_type), intent(inout)  zth,
type(prognostic_field_type), intent(inout)  th,
type(prognostic_field_type), dimension(:), intent(inout)  zq,
type(prognostic_field_type), dimension(:), intent(inout)  q 
)
private

Definition at line 196 of file lowerbc.F90.

196  type(model_state_type), target, intent(inout) :: current_state
197  type(prognostic_field_type), intent(inout) :: zu, zv, th, zth
198  type(prognostic_field_type), dimension(:), intent(inout) :: q, zq
199  integer, intent(in) :: current_y_index, current_x_index
200 
201  integer :: n
202  real(kind=DEFAULT_PRECISION) :: horizontal_velocity_at_k2
203 
204  if (.not. current_state%use_viscosity_and_diffusion .or. .not. current_state%use_surface_boundary_conditions) then
205  call simple_boundary_values(current_state, current_y_index, current_x_index, th, q)
206  else
207  ! one level in for the halo-column
208  if (.not. (current_state%column_local_y .lt. current_state%local_grid%local_domain_start_index(y_index)-1 .or.&
209  current_state%column_local_x .lt. current_state%local_grid%local_domain_start_index(x_index)-1 .or.&
210  current_state%column_local_y .gt. current_state%local_grid%local_domain_end_index(y_index)+1 .or.&
211  current_state%column_local_x .gt. current_state%local_grid%local_domain_end_index(x_index)+1)) then
212 
213  !if (.not. current_state%halo_column) then
214  ! Include one halo to ensure that the halo is set in tvdadvection. This is done using the
215  ! logic from the timestep callback in tvdadvection in the timestep callback above
216  horizontal_velocity_at_k2=0.0_default_precision
217 #ifdef U_ACTIVE
218  horizontal_velocity_at_k2=(0.5_default_precision*(current_state%zu%data(2,current_y_index,current_x_index)+&
219  zu%data(2,current_y_index,current_x_index-1))+ current_state%ugal)**2
220 #endif
221 #ifdef V_ACTIVE
222  horizontal_velocity_at_k2=horizontal_velocity_at_k2+(0.5_default_precision*(zv%data(&
223  2,current_y_index,current_x_index)+zv%data(2,current_y_index-1,current_x_index))+current_state%vgal)**2
224 #endif
225  horizontal_velocity_at_k2=sqrt(horizontal_velocity_at_k2)+smallp
226 
227  if (current_state%type_of_surface_boundary_conditions == prescribed_surface_fluxes) then
228  call compute_using_fixed_surface_fluxes(current_state, current_y_index, current_x_index, &
229  horizontal_velocity_at_k2, th, q)
230  else
231  call compute_using_fixed_surface_temperature(current_state, current_y_index, current_x_index, &
232  horizontal_velocity_at_k2, zth, th, zq, q)
233  end if
234 
235  current_state%dis%data(1, current_y_index, current_x_index)=0.0_default_precision
236  current_state%dis_th%data(1, current_y_index, current_x_index)=0.0_default_precision
237 
238  if (current_state%backscatter) then
239  do n=1, current_state%number_q_fields
240  current_state%disq(n)%data(1,current_y_index,current_x_index)=0.0_default_precision
241  end do
242  end if
243 
244  !-----------------------
245  ! _return viscous number
246  !-----------------------
247 
248  current_state%cvis=max(current_state%cvis,max(current_state%vis_coefficient%data(1, current_y_index, current_x_index),&
249  current_state%diff_coefficient%data(1, current_y_index, current_x_index))*viscous_courant_coefficient)
250  ! CVIS will be multiplied by DTM_X4 in TESTCFL
251  else if (current_x_index == 1 .and. current_y_index == 1) then
252  call register_async_wrapping_recv_requests(current_state)
253  else if (current_x_index == current_state%local_grid%local_domain_end_index(x_index)+&
254  current_state%local_grid%halo_size(x_index) .and. current_y_index == &
255  current_state%local_grid%local_domain_end_index(y_index)+current_state%local_grid%halo_size(y_index)) then
256  call complete_async_wrapping(current_state, zth, zq)
257  end if
258  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ compute_using_fixed_surface_fluxes()

subroutine lowerbc_mod::compute_using_fixed_surface_fluxes ( type(model_state_type), intent(inout), target  current_state,
integer, intent(in)  current_y_index,
integer, intent(in)  current_x_index,
real(kind=default_precision)  horizontal_velocity_at_k2,
type(prognostic_field_type), intent(inout)  th,
type(prognostic_field_type), dimension(:), intent(inout)  q 
)
private

Definition at line 631 of file lowerbc.F90.

631  type(model_state_type), target, intent(inout) :: current_state
632  type(prognostic_field_type), intent(inout) :: th
633  type(prognostic_field_type), dimension(:), intent(inout) :: q
634  integer, intent(in) :: current_y_index, current_x_index
635  real(kind=DEFAULT_PRECISION) :: horizontal_velocity_at_k2
636 
637  if (current_state%fbuoy .gt. 0.0_default_precision) then
638  call handle_convective_fluxes(current_state, current_y_index, current_x_index, horizontal_velocity_at_k2, th, q)
639  else if (current_state%fbuoy .eq. 0.0_default_precision) then
640  call handle_neutral_fluxes(current_state, current_y_index, current_x_index, horizontal_velocity_at_k2, th, q)
641  else
642  call handle_stable_fluxes(current_state, current_y_index, current_x_index, horizontal_velocity_at_k2, th, q)
643  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ compute_using_fixed_surface_temperature()

subroutine lowerbc_mod::compute_using_fixed_surface_temperature ( type(model_state_type), intent(inout), target  current_state,
integer, intent(in)  current_y_index,
integer, intent(in)  current_x_index,
real(kind=default_precision), intent(in)  horizontal_velocity_at_k2,
type(prognostic_field_type), intent(inout)  zth,
type(prognostic_field_type), intent(inout)  th,
type(prognostic_field_type), dimension(:), intent(inout)  zq,
type(prognostic_field_type), dimension(:), intent(inout)  q 
)
private

Definition at line 649 of file lowerbc.F90.

649  type(model_state_type), target, intent(inout) :: current_state
650  type(prognostic_field_type), intent(inout) :: th, zth
651  type(prognostic_field_type), dimension(:), intent(inout) :: q, zq
652  real(kind=DEFAULT_PRECISION), intent(in) :: horizontal_velocity_at_k2
653  integer, intent(in) :: current_y_index, current_x_index
654 
655  real(kind=DEFAULT_PRECISION) :: dthv_surf, ustr, thvstr
656  integer :: convergence_status, n
657 
658  if (current_state%passive_q) then ! i.e. q is not active
659  ! Assuming no liquid water at level 2
660  dthv_surf = zth%data(2, current_y_index, current_x_index) + &
661  current_state%global_grid%configuration%vertical%thref(2) - current_state%theta_virtual_surf
662  else
663  dthv_surf=zth%data(2, current_y_index, current_x_index) + current_state%global_grid%configuration%vertical%thref(2)&
664  *(1.0_default_precision + current_state%cq(current_state%water_vapour_mixing_ratio_index)*&
665  zq(current_state%water_vapour_mixing_ratio_index)%data(2,current_y_index,current_x_index)) - &
666  current_state%theta_virtual_surf
667  end if
668  convergence_status = mostbc(current_state, horizontal_velocity_at_k2, dthv_surf,&
669  current_state%global_grid%configuration%vertical%zn(2), ustr, thvstr)
670 
671  current_state%vis_coefficient%data(1, current_y_index, current_x_index)=&
672  current_state%global_grid%configuration%vertical%czn*ustr**2/horizontal_velocity_at_k2
673  current_state%diff_coefficient%data(1, current_y_index, current_x_index)=&
674  current_state%global_grid%configuration%vertical%czn*ustr*thvstr/(dthv_surf+smallp)
675  zth%data(1, current_y_index, current_x_index)=2.0*current_state%theta_surf-zth%data(2, current_y_index, current_x_index)-&
676  (current_state%global_grid%configuration%vertical%thref(2)+current_state%global_grid%configuration%vertical%thref(1))
677  th%data(1, current_y_index, current_x_index)=2.0*current_state%theta_surf-th%data(2, current_y_index, current_x_index)-&
678  (current_state%global_grid%configuration%vertical%thref(2)+current_state%global_grid%configuration%vertical%thref(1))
679 
680  if (current_state%number_q_fields .gt. 0) then
681  n=iqv
682  zq(n)%data(1, current_y_index, current_x_index)=zq(n)%data(2, current_y_index, current_x_index)
683  q(n)%data(1, current_y_index, current_x_index)=q(n)%data(2, current_y_index, current_x_index)
684  if (.not. current_state%passive_q) then
685  zq(current_state%water_vapour_mixing_ratio_index)%data(1,current_y_index,current_x_index)=&
686  2.0_default_precision*current_state%surface_vapour_mixing_ratio-&
687  zq(current_state%water_vapour_mixing_ratio_index)%data(2,current_y_index,current_x_index)
688  q(current_state%water_vapour_mixing_ratio_index)%data(1,current_y_index,current_x_index)=&
689  2.0_default_precision*current_state%surface_vapour_mixing_ratio-&
690  q(current_state%water_vapour_mixing_ratio_index)%data(2,current_y_index,current_x_index)
691  endif
692  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ finalisation_callback()

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

Definition at line 148 of file lowerbc.F90.

148  type(model_state_type), target, intent(inout) :: current_state
149 
150  if (allocated(x_wrapping_send_buffer)) deallocate(x_wrapping_send_buffer)
151  if (allocated(y_wrapping_send_buffer)) deallocate(y_wrapping_send_buffer)
152  if (allocated(x_wrapping_recv_buffer)) deallocate(x_wrapping_recv_buffer)
153  if (allocated(y_wrapping_recv_buffer)) deallocate(y_wrapping_recv_buffer)
Here is the caller graph for this function:

◆ handle_convective_fluxes()

subroutine lowerbc_mod::handle_convective_fluxes ( type(model_state_type), intent(inout), target  current_state,
integer, intent(in)  current_y_index,
integer, intent(in)  current_x_index,
real(kind=default_precision), intent(in)  horizontal_velocity_at_k2,
type(prognostic_field_type), intent(inout)  th,
type(prognostic_field_type), dimension(:), intent(inout)  q 
)
private

Definition at line 493 of file lowerbc.F90.

493  type(model_state_type), target, intent(inout) :: current_state
494  type(prognostic_field_type), intent(inout) :: th
495  type(prognostic_field_type), dimension(:), intent(inout) :: q
496  integer, intent(in) :: current_y_index, current_x_index
497  real(kind=DEFAULT_PRECISION), intent(in) :: horizontal_velocity_at_k2
498 
499  integer :: n
500  real(kind=DEFAULT_PRECISION) :: ustr
501 
502  ustr=look(current_state, horizontal_velocity_at_k2)
503 
504  current_state%vis_coefficient%data(1, current_y_index, current_x_index)=&
505  current_state%global_grid%configuration%vertical%czn*ustr**2/ horizontal_velocity_at_k2
506  current_state%diff_coefficient%data(1, current_y_index, current_x_index)=&
507  (von_karman_constant*current_state%global_grid%configuration%vertical%czn*ustr/alphah)/&
508  (current_state%global_grid%configuration%vertical%zlogth- 2.*log(&
509  (1.+sqrt(1.+gammah*von_karman_constant*current_state%fbuoy*(current_state%global_grid%configuration%vertical%zn(2)+z0)&
510  /ustr**3))/ (1.+sqrt(1.+gammah*von_karman_constant*current_state%fbuoy*z0th/ustr**3))))
511  if (current_state%th%active) th%data(1, current_y_index, current_x_index)=th%data(2, current_y_index, current_x_index)+&
512  current_state%surface_temperature_flux*current_state%global_grid%configuration%vertical%dzn(2)/&
513  current_state%diff_coefficient%data(1, current_y_index, current_x_index)
514 
515 
516  ! Surface Flux of vapour
517  if (current_state%number_q_fields .gt. 0) then
518  n=iqv
519  q(n)%data(1, current_y_index, current_x_index)=q(n)%data(2, current_y_index, current_x_index)+&
520  current_state%surface_vapour_flux*current_state%global_grid%configuration%vertical%dzn(2)/&
521  current_state%diff_coefficient%data(1, current_y_index, current_x_index)
522  endif
Here is the call graph for this function:
Here is the caller graph for this function:

◆ handle_neutral_fluxes()

subroutine lowerbc_mod::handle_neutral_fluxes ( type(model_state_type), intent(inout), target  current_state,
integer, intent(in)  current_y_index,
integer, intent(in)  current_x_index,
real(kind=default_precision)  horizontal_velocity_at_k2,
type(prognostic_field_type), intent(inout)  th,
type(prognostic_field_type), dimension(:), intent(inout)  q 
)
private

Definition at line 551 of file lowerbc.F90.

551  type(model_state_type), target, intent(inout) :: current_state
552  type(prognostic_field_type), intent(inout) :: th
553  type(prognostic_field_type), dimension(:), intent(inout) :: q
554  integer, intent(in) :: current_y_index, current_x_index
555  real(kind=DEFAULT_PRECISION) :: horizontal_velocity_at_k2
556 
557  integer :: n
558  real(kind=DEFAULT_PRECISION) :: ustr
559 
560  ustr=horizontal_velocity_at_k2*current_state%global_grid%configuration%vertical%vk_on_zlogm
561  current_state%vis_coefficient%data(1, current_y_index, current_x_index)=current_state%global_grid%configuration%vertical%czn*&
562  ustr**2/horizontal_velocity_at_k2
563  current_state%diff_coefficient%data(1, current_y_index, current_x_index)=&
564  current_state%vis_coefficient%data(1, current_y_index, current_x_index)*&
565  current_state%global_grid%configuration%vertical%zlogm/(alphah*current_state%global_grid%configuration%vertical%zlogth)
566  if (current_state%th%active) th%data(1, current_y_index, current_x_index)=th%data(2, current_y_index, current_x_index)+&
567  current_state%surface_temperature_flux*current_state%global_grid%configuration%vertical%dzn(2)/&
568  current_state%diff_coefficient%data(1, current_y_index, current_x_index)
569 
570  ! Flux of vapour
571  if (current_state%number_q_fields .gt. 0) then
572  n=iqv
573  q(n)%data(1, current_y_index, current_x_index)=q(n)%data(2, current_y_index, current_x_index)+&
574  current_state%surface_vapour_flux*current_state%global_grid%configuration%vertical%dzn(2)/&
575  current_state%diff_coefficient%data(1, current_y_index, current_x_index)
576  endif
Here is the caller graph for this function:

◆ handle_stable_fluxes()

subroutine lowerbc_mod::handle_stable_fluxes ( type(model_state_type), intent(inout), target  current_state,
integer, intent(in)  current_y_index,
integer, intent(in)  current_x_index,
real(kind=default_precision)  horizontal_velocity_at_k2,
type(prognostic_field_type), intent(inout)  th,
type(prognostic_field_type), dimension(:), intent(inout)  q 
)
private

Definition at line 580 of file lowerbc.F90.

580  type(model_state_type), target, intent(inout) :: current_state
581  type(prognostic_field_type), intent(inout) :: th
582  type(prognostic_field_type), dimension(:), intent(inout) :: q
583  integer, intent(in) :: current_y_index, current_x_index
584  real(kind=DEFAULT_PRECISION) :: horizontal_velocity_at_k2
585 
586  real(kind=DEFAULT_PRECISION) :: ustr
587  integer :: n
588 
589  ustr=horizontal_velocity_at_k2*current_state%global_grid%configuration%vertical%vk_on_zlogm
590  if((current_state%fbuoy - 1.e-9_default_precision) .lt. -4.0_default_precision*&
591  von_karman_constant**2*horizontal_velocity_at_k2**3/ (27.0_default_precision*betam*&
592  current_state%global_grid%configuration%vertical%zn(2)*current_state%global_grid%configuration%vertical%zlogm**2)) then
593  ! Too stable for turbulence
594  current_state%vis_coefficient%data(1, current_y_index, current_x_index)=0.0_default_precision
595  current_state%diff_coefficient%data(1, current_y_index, current_x_index)=0.0_default_precision
596  ! Level 1 values of l_th and l_q set to be harmless for advection scheme
597  if (current_state%th%active) th%data(1, current_y_index, current_x_index)=th%data(2, current_y_index, current_x_index)
598  if (current_state%number_q_fields .gt. 0) then
599  do n=1, current_state%number_q_fields
600  q(n)%data(1, current_y_index, current_x_index)=q(n)%data(2, current_y_index, current_x_index)
601  end do
602  end if
603  else
604  ustr=ustr/3.0_default_precision*(1.0_default_precision-2.0_default_precision*&
605  cos((acos(-27.0_default_precision*betam*von_karman_constant*current_state%global_grid%configuration%vertical%zn(2)*&
606  current_state%fbuoy/(current_state%global_grid%configuration%vertical%zlogm*&
607  2.0_default_precision*ustr**3)-1.0_default_precision)+ 2.0_default_precision*pi)/3.0_default_precision))
608  current_state%vis_coefficient%data(1, current_y_index, current_x_index)=&
609  current_state%global_grid%configuration%vertical%czn*ustr**2/horizontal_velocity_at_k2
610  current_state%diff_coefficient%data(1, current_y_index, current_x_index)=&
611  current_state%vis_coefficient%data(1, current_y_index, current_x_index)*&
612  (current_state%global_grid%configuration%vertical%zlogm-betam*current_state%global_grid%configuration%vertical%zn(2)*&
613  von_karman_constant*current_state%fbuoy/ustr**3)/(alphah*current_state%global_grid%configuration%vertical%zlogth-betah*&
614  von_karman_constant*current_state%fbuoy* (current_state%global_grid%configuration%vertical%zn(2)+ z0-z0th)/ustr**3)
615  if (current_state%th%active) th%data(1, current_y_index, current_x_index)=th%data(2, current_y_index, current_x_index)+&
616  current_state%surface_temperature_flux*current_state%global_grid%configuration%vertical%dzn(2)/&
617  current_state%diff_coefficient%data(1, current_y_index, current_x_index)
618 
619  !Flux of vapour
620  if (current_state%number_q_fields .gt. 0) then
621  n=iqv
622  q(n)%data(1, current_y_index, current_x_index)=q(n)%data(2, current_y_index, current_x_index)+&
623  current_state%surface_vapour_flux*current_state%global_grid%configuration%vertical%dzn(2)/&
624  current_state%diff_coefficient%data(1, current_y_index, current_x_index)
625  endif
626  end if
Here is the caller graph for this function:

◆ initialisation_callback()

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

Definition at line 50 of file lowerbc.F90.

50  type(model_state_type), target, intent(inout) :: current_state
51 
52  real(kind=DEFAULT_PRECISION) :: bhbc
53  integer :: num_wrapped_fields
54 
55  ! Adhill - this check is only required so that the vis_ and diff_coefficients
56  ! are allocated in their respective components
57  if (.not. is_component_enabled(current_state%options_database, "diffusion")) then
58  call log_master_log(log_error, "Lowerbc requires the diffusion component to be enabled")
59  end if
60  if (.not. is_component_enabled(current_state%options_database, "viscosity")) then
61  call log_master_log(log_error, "Lowerbc requires the viscosity component to be enabled")
62  end if
63 
64  call allocate_applicable_fields(current_state)
65 
66  wrapping_comm_requests=mpi_request_null
67 
68  num_wrapped_fields=0
69  if (current_state%th%active) num_wrapped_fields=1
70  num_wrapped_fields=num_wrapped_fields+current_state%number_q_fields
71 
72  if (num_wrapped_fields .gt. 0) then
73  if (current_state%parallel%my_coords(y_index) == 0 .or. &
74  current_state%parallel%my_coords(y_index) == current_state%parallel%dim_sizes(y_index)-1) then
75  if (current_state%parallel%my_coords(y_index) == 0) then
76  y_wrapping_target_id=current_state%local_grid%neighbours(y_index, 1)
77  else
78  y_wrapping_target_id=current_state%local_grid%neighbours(y_index, 3)
79  end if
80  if (current_state%parallel%my_rank .ne. y_wrapping_target_id) then
81  allocate(y_wrapping_send_buffer(current_state%local_grid%size(x_index)+4, 2, num_wrapped_fields), &
82  y_wrapping_recv_buffer(current_state%local_grid%size(x_index)+4, 2, num_wrapped_fields))
83  end if
84  end if
85 
86  if (current_state%parallel%my_coords(x_index) == 0 .or. &
87  current_state%parallel%my_coords(x_index) == current_state%parallel%dim_sizes(x_index)-1) then
88  if (current_state%parallel%my_coords(x_index) == 0) then
89  x_wrapping_target_id=current_state%local_grid%neighbours(x_index, 1)
90  else
91  x_wrapping_target_id=current_state%local_grid%neighbours(x_index, 3)
92  end if
93  if (current_state%parallel%my_rank .ne. x_wrapping_target_id) then
94  allocate(x_wrapping_send_buffer(current_state%local_grid%size(y_index)+4, 2, num_wrapped_fields), &
95  x_wrapping_recv_buffer(current_state%local_grid%size(y_index)+4, 2, num_wrapped_fields))
96  end if
97  end if
98  end if
99 
100  viscous_courant_coefficient=1.0_default_precision/current_state%global_grid%configuration%vertical%dzn(2)**2+&
101  1.0_default_precision/(current_state%global_grid%configuration%horizontal%dx*&
102  current_state%global_grid%configuration%horizontal%dx)+&
103  1.0_default_precision/(current_state%global_grid%configuration%horizontal%dy*&
104  current_state%global_grid%configuration%horizontal%dy)
105 
106  if ( current_state%use_surface_boundary_conditions .and. &
107  current_state%type_of_surface_boundary_conditions == prescribed_surface_values) then
108  ! variables below are only required when PRESCRIBED_SURFACE_VALUES are used.
109  tstrcona=von_karman_constant/alphah*current_state%global_grid%configuration%vertical%zlogth
110  bhbc=alphah*current_state%global_grid%configuration%vertical%zlogth
111  rhmbc=betah*(current_state%global_grid%configuration%vertical%zn(2)+z0-z0th)/&
112  (betam*current_state%global_grid%configuration%vertical%zn(2))
113  ddbc=current_state%global_grid%configuration%vertical%zlogm*(bhbc-&
114  rhmbc*current_state%global_grid%configuration%vertical%zlogm)
115  ddbc_x4=4.*ddbc
116  r2ddbc=0.5_default_precision/ddbc
117  eecon=2.0_default_precision*rhmbc*current_state%global_grid%configuration%vertical%zlogm-bhbc
118  rcmbc=1.0_default_precision/current_state%cmbc
119  tstrconb=von_karman_constant/alphah
120  x4con=gammam*(current_state%global_grid%configuration%vertical%zn(2)+z0)
121  xx0con=gammam*z0
122  y2con=gammah*(current_state%global_grid%configuration%vertical%zn(2)+z0)
123  yy0con=gammah*z0th
124  else
125  tstrcona=0.0
126  bhbc=0.0
127  rhmbc=0.0
128  ddbc=0.0
129  ddbc_x4=0.0
130  r2ddbc=0.0
131  eecon=0.0
132  rcmbc=0.0
133  tstrconb=0.0
134  x4con=0.0
135  xx0con=0.0
136  y2con=0.0
137  yy0con=0.0
138  endif
139 
140  ! Determine vapour index
141  if (.not. current_state%passive_q) then
142  iqv = get_q_index(standard_q_names%VAPOUR, 'lowerbc')
143  endif
144 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ look()

real(kind=default_precision) function lowerbc_mod::look ( type(model_state_type), intent(inout), target  current_state,
real(kind=default_precision), intent(in)  vel 
)
private

Definition at line 526 of file lowerbc.F90.

526  type(model_state_type), target, intent(inout) :: current_state
527  real(kind=DEFAULT_PRECISION), intent(in) :: vel ! Horizontal speed at lowest model level
528 
529  real(kind=DEFAULT_PRECISION) :: lookup_real_posn
530  integer :: lookup_int_posn
531 
532  lookup_real_posn=1.0_default_precision+real(current_state%lookup_table_entries-1)*&
533  log(vel/current_state%velmin)*current_state%aloginv
534  lookup_int_posn=int(lookup_real_posn)
535 
536  if (lookup_int_posn .ge. 1) then
537  if (lookup_int_posn .lt. current_state%lookup_table_entries) then ! Linear interpolation
538  look=current_state%lookup_table_ustr(lookup_int_posn)+ (lookup_real_posn-real(lookup_int_posn))*&
539  (current_state%lookup_table_ustr(lookup_int_posn+1)-current_state%lookup_table_ustr(lookup_int_posn))
540  else ! Near neutral
541  look=vel*current_state%cneut
542  end if
543  else ! Nearly free convection
544  look=vel**(-convective_limit)*current_state%cfc
545  end if
Here is the caller graph for this function:

◆ lowerbc_get_descriptor()

type(component_descriptor_type) function, public lowerbc_mod::lowerbc_get_descriptor ( )

Descriptor of this component for registration.

Returns
The diverr component descriptor

Definition at line 42 of file lowerbc.F90.

42  lowerbc_get_descriptor%name="lower_bc"
43  lowerbc_get_descriptor%version=0.1
44  lowerbc_get_descriptor%initialisation=>initialisation_callback
45  lowerbc_get_descriptor%timestep=>timestep_callback
46  lowerbc_get_descriptor%finalisation=>finalisation_callback
Here is the call graph for this function:

◆ mostbc()

integer function lowerbc_mod::mostbc ( type(model_state_type), intent(inout), target  current_state,
real(kind=default_precision), intent(in)  delu,
real(kind=default_precision), intent(in)  delt,
real(kind=default_precision), intent(in)  z1,
real(kind=default_precision), intent(out)  ustrdg,
real(kind=default_precision), intent(out)  tstrdg 
)
private

Solves the Monin-Obukhov equations in the case of specified surface values of temperature and mixing ratio,combined into a specified value of virtual temperature. It is a modified version of the subroutine described in Bull and Derbyshire (TDN197) based on the assumption that the similarity functions and roughness lengths for temperature and mixing ratio are the same. In that case, all the original theory can be used if we replace temperature by virtual temperature. The form of the non-dimensionalised wind shear used is 1.0 + BETAM z/L for the stable case and (1.0 - GAMMAM z/L)**1/4 for the unstable case. The form of the non-dimensionalised temperature gradient used is 1.0 + BETAH z/L for the stable case and (1.0 - GAMMAH z/L)**1/2 for the unstable case.

Parameters
current_stateThe current model state
deluThe wind speed at the lowest grid point
deltThe virtual potential temperature difference between the lowest grid point and the surface
z1The height of the lowest grid point ABOVE the roughness length Z0
ustrdgThe diagnosed value of friction velocity
tstrdgThe diagnosed value of surface virtual temperature scale
Returns
The convergence criteria - success, richardson number too large or failure

Definition at line 735 of file lowerbc.F90.

735  type(model_state_type), target, intent(inout) :: current_state
736  real(kind=DEFAULT_PRECISION), intent(in) :: delu, delt, z1
737  real(kind=DEFAULT_PRECISION), intent(out) :: ustrdg, tstrdg
738 
739  if (delt .lt. 0.0_default_precision) then
740  if (delu .le. smallp) then
741  ustrdg=0.0_default_precision
742  tstrdg=tstrcona*delt
743  else
744  ! The unstable case
745  mostbc=solve_monin_obukhov_unstable_case(delu, delt, current_state%ellmocon, &
746  ustrdg, tstrdg, current_state%global_grid%configuration%vertical)
747  end if
748  else if (delt .gt. 0.0_default_precision) then
749  ! The stable case
750  mostbc=solve_monin_obukhov_stable_case(delu, delt, current_state%global_grid%configuration%vertical%zlogm, &
751  current_state%cmbc, ustrdg, tstrdg)
752  else
753  ! Trivial neutral case
754  ustrdg=current_state%global_grid%configuration%vertical%vk_on_zlogm*delu
755  tstrdg=0.0_default_precision
756  mostbc=convergence_success
757  end if
758 
759  if (mostbc .ne. convergence_success) then
760  if (mostbc .eq. convergence_richardson_too_large) then
761  call log_log(log_warn, "Richardson number greater than critical value")
762  else if(mostbc .eq. convergence_failure) then
763  call log_log(log_error, "Convergence failure after 200 iterations")
764  end if
765  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ package_x_wrapping_send_buffer()

subroutine lowerbc_mod::package_x_wrapping_send_buffer ( type(model_state_type), intent(inout), target  current_state,
type(prognostic_field_type), intent(inout)  zth,
type(prognostic_field_type), dimension(:), intent(inout)  zq,
integer, intent(in)  first_x_index,
integer, intent(in)  second_x_index 
)
private

Packages theta and Q fields (if enabled) into the send buffer for X.

Parameters
current_stateThe current model state
zthTemperature field
zqQ fields
first_x_indexThe first X index to read from the data field
second_x_indexThe second X index to read from the data field

Definition at line 415 of file lowerbc.F90.

415  type(model_state_type), target, intent(inout) :: current_state
416  type(prognostic_field_type), intent(inout) :: zth
417  type(prognostic_field_type), dimension(:), intent(inout) :: zq
418  integer, intent(in) :: first_x_index, second_x_index
419 
420  integer :: index_start, n
421 
422  index_start=0
423  if (current_state%th%active) then
424  x_wrapping_send_buffer(:,1,1)=zth%data(1,:,first_x_index)
425  x_wrapping_send_buffer(:,2,1)=zth%data(1,:,second_x_index)
426  index_start=index_start+1
427  end if
428  if (current_state%number_q_fields .gt. 0) then
429  do n=1, current_state%number_q_fields
430  x_wrapping_send_buffer(:,1,index_start+n)= zq(n)%data(1,:,first_x_index)
431  x_wrapping_send_buffer(:,2,index_start+n)= zq(n)%data(1,:,second_x_index)
432  end do
433  end if
Here is the caller graph for this function:

◆ package_y_wrapping_send_buffer()

subroutine lowerbc_mod::package_y_wrapping_send_buffer ( type(model_state_type), intent(inout), target  current_state,
type(prognostic_field_type), intent(inout)  zth,
type(prognostic_field_type), dimension(:), intent(inout)  zq,
integer, intent(in)  first_y_index,
integer, intent(in)  second_y_index 
)
private

Packages theta and Q fields (if enabled) into the send buffer for Y.

Parameters
current_stateThe current model state
zthTemperature field
zqQ fields
first_y_indexThe first Y index to read from the data field
second_y_indexThe second Y index to read from the data field

Definition at line 387 of file lowerbc.F90.

387  type(model_state_type), target, intent(inout) :: current_state
388  type(prognostic_field_type), intent(inout) :: zth
389  type(prognostic_field_type), dimension(:), intent(inout) :: zq
390  integer, intent(in) :: first_y_index, second_y_index
391 
392  integer :: index_start, n
393 
394  index_start=0
395  if (current_state%th%active) then
396  y_wrapping_send_buffer(:,1,1)=zth%data(1,first_y_index,:)
397  y_wrapping_send_buffer(:,2,1)=zth%data(1,second_y_index,:)
398  index_start=index_start+1
399  end if
400  if (current_state%number_q_fields .gt. 0) then
401  do n=1, current_state%number_q_fields
402  y_wrapping_send_buffer(:,1,index_start+n)=zq(n)%data(1,first_y_index,:)
403  y_wrapping_send_buffer(:,2,index_start+n)=zq(n)%data(1,second_y_index,:)
404  end do
405  end if
Here is the caller graph for this function:

◆ register_async_wrapping_recv_requests()

subroutine lowerbc_mod::register_async_wrapping_recv_requests ( type(model_state_type), intent(inout), target  current_state)
private

Registers asynchronous wrapping recv requests as needed.

Parameters
current_stateThe current model state

Definition at line 264 of file lowerbc.F90.

264  type(model_state_type), target, intent(inout) :: current_state
265 
266  integer :: ierr
267 
268  if (allocated(y_wrapping_recv_buffer)) then
269  call mpi_irecv(y_wrapping_recv_buffer, size(y_wrapping_recv_buffer), precision_type, &
270  y_wrapping_target_id, 0, current_state%parallel%neighbour_comm, wrapping_comm_requests(1), ierr)
271  end if
272  if (allocated(x_wrapping_recv_buffer)) then
273  call mpi_irecv(x_wrapping_recv_buffer, size(x_wrapping_recv_buffer), precision_type, &
274  x_wrapping_target_id, 0, current_state%parallel%neighbour_comm, wrapping_comm_requests(3), ierr)
275  end if
Here is the caller graph for this function:

◆ simple_boundary_values()

subroutine lowerbc_mod::simple_boundary_values ( type(model_state_type), intent(inout), target  current_state,
integer, intent(in)  current_y_index,
integer, intent(in)  current_x_index,
type(prognostic_field_type), intent(inout)  th,
type(prognostic_field_type), dimension(:), intent(inout)  q 
)
private

Definition at line 696 of file lowerbc.F90.

696  type(model_state_type), target, intent(inout) :: current_state
697  type(prognostic_field_type), intent(inout) :: th
698  type(prognostic_field_type), dimension(:), intent(inout) :: q
699  integer, intent(in) :: current_y_index, current_x_index
700 
701  integer :: n
702 
703  current_state%vis_coefficient%data(1, current_y_index, current_x_index)=0.0_default_precision
704  current_state%diff_coefficient%data(1, current_y_index, current_x_index)=0.0_default_precision
705  if(current_state%backscatter) current_state%dis%data(1, current_y_index, current_x_index)=0.0_default_precision
706  if(current_state%backscatter .and. current_state%th%active) &
707  current_state%dis_th%data(1, current_y_index, current_x_index)=0.0_default_precision
708  if (current_state%th%active) then
709  th%data(1, current_y_index, current_x_index) = th%data(2, current_y_index, current_x_index)
710  end if
711  if (current_state%number_q_fields .gt. 0) then
712  do n=1, current_state%number_q_fields
713  q(n)%data(1, current_y_index, current_x_index) = q(n)%data(2, current_y_index, current_x_index)
714  if (current_state%backscatter) current_state%disq(n)%data(1, current_y_index, current_x_index) = 0.0_default_precision
715  end do
716  end if
Here is the caller graph for this function:

◆ solve_monin_obukhov_stable_case()

integer function lowerbc_mod::solve_monin_obukhov_stable_case ( real(kind=default_precision), intent(in)  delu,
real(kind=default_precision), intent(in)  delt,
real(kind=default_precision), intent(in)  zlogm,
real(kind=default_precision), intent(in)  cmbc,
real(kind=default_precision), intent(out)  ustrdg,
real(kind=default_precision), intent(out)  tstrdg 
)
private

Definition at line 819 of file lowerbc.F90.

819  real(kind=DEFAULT_PRECISION), intent(in) :: delu, delt, zlogm, cmbc
820  real(kind=DEFAULT_PRECISION), intent(out) :: ustrdg, tstrdg
821 
822  real(kind=DEFAULT_PRECISION) :: am, ah, ee, ff, det
823 
824  am=von_karman_constant*delu
825  ah=von_karman_constant*delt
826  ee=am*eecon
827  ff=ah*cmbc-rhmbc*am*am !
828  det=ee*ee-ddbc_x4*ff
829  solve_monin_obukhov_stable_case=convergence_success
830  ! Test for laminar flow
831  if (ff .gt. 0.0_default_precision) then
832  if ((ee .lt. 0.0_default_precision).and.(det .gt. 0.0_default_precision)) then
833  ustrdg=(-ee+sqrt(det))*r2ddbc
834  tstrdg=ustrdg*(am-zlogm*ustrdg)*rcmbc
835  else
836  solve_monin_obukhov_stable_case=convergence_richardson_too_large
837  ustrdg=0.0_default_precision
838  tstrdg=0.0_default_precision
839  end if
840  else if (ddbc .eq. 0.0_default_precision) then
841  ! Degenerate case
842  ustrdg=-ff/ee
843  tstrdg=delt*ustrdg/delu
844  else
845  ! Solve quadratic for USTRDG
846  ustrdg=(-ee+sqrt(det))*r2ddbc
847  tstrdg=ustrdg*(am-zlogm*ustrdg)*rcmbc
848  end if
Here is the caller graph for this function:

◆ solve_monin_obukhov_unstable_case()

integer function lowerbc_mod::solve_monin_obukhov_unstable_case ( real(kind=default_precision), intent(in)  delu,
real(kind=default_precision), intent(in)  delt,
real(kind=default_precision), intent(in)  ellmocon,
real(kind=default_precision), intent(out)  ustrdg,
real(kind=default_precision), intent(out)  tstrdg,
type(vertical_grid_configuration_type), intent(inout)  vertical_grid 
)
private

Definition at line 769 of file lowerbc.F90.

769  real(kind=DEFAULT_PRECISION), intent(in) :: delu, delt, ellmocon
770  real(kind=DEFAULT_PRECISION), intent(out) :: ustrdg, tstrdg
771  type(vertical_grid_configuration_type), intent(inout) :: vertical_grid
772 
773  integer :: i
774  real(kind=DEFAULT_PRECISION) :: ellmo, psim, psih, x4, xx, xx0, y2, yy, yy0, err_ustr, err_tstr, &
775  ustrl, tstrl, & ! U and T star at start of iteration
776  ustrpl, tstrpl ! U and T star at end of iteration
777 
778  ! First set initial values
779  ustrl=vertical_grid%vk_on_zlogm*delu
780  tstrl=tstrcona*delt
781 
782  ! Now start iteration
783  do i=1, 200
784  ellmo=ustrl*ustrl*ellmocon/tstrl
785 
786  ! Test for possible square root of negative quantity
787  x4=1.0_default_precision-x4con/ellmo
788  if (x4 .lt. 0.0_default_precision) call log_log(log_error, "Negative square root in x4")
789 
790  xx=sqrt(sqrt(x4))
791  xx0=sqrt(sqrt(1.0_default_precision-xx0con / ellmo))
792  psim=2.*( log((xx+1.0_default_precision)/(xx0+1.0_default_precision))-atan(xx)+atan(xx0) )+&
793  log((xx*xx+1.0_default_precision)/(xx0*xx0+1.0_default_precision))
794  ustrpl=von_karman_constant*delu/(vertical_grid%zlogm-psim)
795 
796  ! Test for possible square root of negative quantity
797  y2=1.-y2con/ellmo
798  if (y2 .lt. 0.0_default_precision) call log_log(log_error, "Negative square root in y2")
799  yy=sqrt(y2)
800  yy0=sqrt(1.0_default_precision-yy0con/ellmo)
801  psih=2.*log((1.0_default_precision+yy)/(1.0_default_precision+yy0))
802  tstrpl=tstrconb*delt/(vertical_grid%zlogth-psih)
803  err_ustr=abs((ustrpl-ustrl)/ ustrl)
804  err_tstr=abs((tstrpl-tstrl)/ tstrl)
805  if ((err_tstr .lt. tolt) .and. (err_ustr .lt. tolm)) then
806  ustrdg=ustrpl
807  tstrdg=tstrpl
808  solve_monin_obukhov_unstable_case=convergence_success
809  return
810  else
811  ustrl=(1.0_default_precision-smth)*ustrpl+smth*ustrl
812  tstrl=(1.0_default_precision-smth)*tstrpl+smth*tstrl
813  end if
814  end do
815  solve_monin_obukhov_unstable_case=convergence_failure
Here is the caller graph for this function:

◆ timestep_callback()

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

Definition at line 174 of file lowerbc.F90.

174  type(model_state_type), target, intent(inout) :: current_state
175 
176  integer :: current_y_index, current_x_index
177 
178  current_y_index=current_state%column_local_y
179  current_x_index=current_state%column_local_x
180 
181  if (current_state%field_stepping == forward_stepping) then
182  call compute_lower_boundary_conditions(current_state, current_y_index, current_x_index, &
183  current_state%u, current_state%v, current_state%th, current_state%th, current_state%q, current_state%q)
184  else
185  if (current_state%scalar_stepping == forward_stepping) then
186  call compute_lower_boundary_conditions(current_state, current_y_index, current_x_index, &
187  current_state%zu, current_state%zv, current_state%th, current_state%zth, current_state%q, current_state%zq)
188  else
189  call compute_lower_boundary_conditions(current_state, current_y_index, current_x_index, &
190  current_state%zu, current_state%zv, current_state%zth, current_state%zth, current_state%zq, current_state%zq)
191  end if
192  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ unpackage_x_wrapping_recv_buffer()

subroutine lowerbc_mod::unpackage_x_wrapping_recv_buffer ( type(model_state_type), intent(inout), target  current_state,
type(prognostic_field_type), intent(inout)  zth,
type(prognostic_field_type), dimension(:), intent(inout)  zq,
integer, intent(in)  first_x_index,
integer, intent(in)  second_x_index 
)
private

Unpackages theta and Q fields from the receive buffer into the fields themselves (if enabled) for X.

Parameters
current_stateThe current model state
zthTemperature field
zqQ fields
first_x_indexThe first X index to read from the data field
second_x_indexThe second X index to read from the data field

Definition at line 471 of file lowerbc.F90.

471  type(model_state_type), target, intent(inout) :: current_state
472  type(prognostic_field_type), intent(inout) :: zth
473  type(prognostic_field_type), dimension(:), intent(inout) :: zq
474  integer, intent(in) :: first_x_index, second_x_index
475 
476  integer :: index_start, n
477 
478  index_start=0
479  if (current_state%th%active) then
480  zth%data(1,:,first_x_index)=x_wrapping_recv_buffer(:,1,1)
481  zth%data(1,:,second_x_index)=x_wrapping_recv_buffer(:,2,1)
482  index_start=index_start+1
483  end if
484  if (current_state%number_q_fields .gt. 0) then
485  do n=1, current_state%number_q_fields
486  zq(n)%data(1,:,first_x_index)=x_wrapping_recv_buffer(:,1,index_start+n)
487  zq(n)%data(1,:,second_x_index)=x_wrapping_recv_buffer(:,2,index_start+n)
488  end do
489  end if
Here is the caller graph for this function:

◆ unpackage_y_wrapping_recv_buffer()

subroutine lowerbc_mod::unpackage_y_wrapping_recv_buffer ( type(model_state_type), intent(inout), target  current_state,
type(prognostic_field_type), intent(inout)  zth,
type(prognostic_field_type), dimension(:), intent(inout)  zq,
integer, intent(in)  first_y_index,
integer, intent(in)  second_y_index 
)
private

Unpackages theta and Q fields from the receive buffer into the fields themselves (if enabled) for Y.

Parameters
current_stateThe current model state
zthTemperature field
zqQ fields
first_y_indexThe first Y index to read from the data field
second_y_indexThe second Y index to read from the data field

Definition at line 443 of file lowerbc.F90.

443  type(model_state_type), target, intent(inout) :: current_state
444  type(prognostic_field_type), intent(inout) :: zth
445  type(prognostic_field_type), dimension(:), intent(inout) :: zq
446  integer, intent(in) :: first_y_index, second_y_index
447 
448  integer :: index_start, n
449 
450  index_start=0
451  if (current_state%th%active) then
452  zth%data(1,first_y_index,:)=y_wrapping_recv_buffer(:,1,1)
453  zth%data(1,second_y_index,:)=y_wrapping_recv_buffer(:,2,1)
454  index_start=index_start+1
455  end if
456  if (current_state%number_q_fields .gt. 0) then
457  do n=1, current_state%number_q_fields
458  zq(n)%data(1,first_y_index,:)=y_wrapping_recv_buffer(:,1,index_start+n)
459  zq(n)%data(1,second_y_index,:)=y_wrapping_recv_buffer(:,2,index_start+n)
460  end do
461  end if
Here is the caller graph for this function:

Variable Documentation

◆ convergence_failure

integer, parameter lowerbc_mod::convergence_failure =3
private

Definition at line 22 of file lowerbc.F90.

◆ convergence_richardson_too_large

integer, parameter lowerbc_mod::convergence_richardson_too_large =2
private

Definition at line 22 of file lowerbc.F90.

◆ convergence_success

integer, parameter lowerbc_mod::convergence_success =1
private

Definition at line 22 of file lowerbc.F90.

22  integer, parameter :: convergence_success=1, convergence_richardson_too_large=2, convergence_failure=3

◆ ddbc

real(kind=default_precision) lowerbc_mod::ddbc
private

Definition at line 27 of file lowerbc.F90.

◆ ddbc_x4

real(kind=default_precision) lowerbc_mod::ddbc_x4
private

Definition at line 27 of file lowerbc.F90.

◆ eecon

real(kind=default_precision) lowerbc_mod::eecon
private

Definition at line 27 of file lowerbc.F90.

◆ iqv

integer lowerbc_mod::iqv
private

Definition at line 33 of file lowerbc.F90.

33  integer :: iqv ! index for vapour

◆ r2ddbc

real(kind=default_precision) lowerbc_mod::r2ddbc
private

Definition at line 27 of file lowerbc.F90.

◆ rcmbc

real(kind=default_precision) lowerbc_mod::rcmbc
private

Definition at line 27 of file lowerbc.F90.

◆ rhmbc

real(kind=default_precision) lowerbc_mod::rhmbc
private

Definition at line 27 of file lowerbc.F90.

◆ smth

real(kind=default_precision), parameter lowerbc_mod::smth = 0.05_DEFAULT_PRECISION
private

Definition at line 24 of file lowerbc.F90.

24  real(kind=DEFAULT_PRECISION), parameter :: smth = 0.05_default_precision,& ! Smoothing between iterations
25  tolm=1.0e-4_default_precision, tolt=1.0e-4_default_precision ! Convergence tollerance for u and t star

◆ tolm

real(kind=default_precision), parameter lowerbc_mod::tolm =1.0E-4_DEFAULT_PRECISION
private

Definition at line 24 of file lowerbc.F90.

◆ tolt

real(kind=default_precision), parameter lowerbc_mod::tolt =1.0E-4_DEFAULT_PRECISION
private

Definition at line 24 of file lowerbc.F90.

◆ tstrcona

real(kind=default_precision) lowerbc_mod::tstrcona
private

Definition at line 27 of file lowerbc.F90.

27  real(kind=DEFAULT_PRECISION) :: tstrcona, rhmbc, ddbc, ddbc_x4, eecon, r2ddbc, rcmbc, tstrconb, &
28  x4con, xx0con, y2con, yy0con, viscous_courant_coefficient

◆ tstrconb

real(kind=default_precision) lowerbc_mod::tstrconb
private

Definition at line 27 of file lowerbc.F90.

◆ viscous_courant_coefficient

real(kind=default_precision) lowerbc_mod::viscous_courant_coefficient
private

Definition at line 27 of file lowerbc.F90.

◆ wrapping_comm_requests

integer, dimension(4) lowerbc_mod::wrapping_comm_requests
private

Definition at line 34 of file lowerbc.F90.

34  integer :: wrapping_comm_requests(4), y_wrapping_target_id, x_wrapping_target_id

◆ x4con

real(kind=default_precision) lowerbc_mod::x4con
private

Definition at line 27 of file lowerbc.F90.

◆ x_wrapping_recv_buffer

real(kind=default_precision), dimension(:,:,:), allocatable lowerbc_mod::x_wrapping_recv_buffer
private

Definition at line 30 of file lowerbc.F90.

◆ x_wrapping_send_buffer

real(kind=default_precision), dimension(:,:,:), allocatable lowerbc_mod::x_wrapping_send_buffer
private

Definition at line 30 of file lowerbc.F90.

30  real(kind=DEFAULT_PRECISION), dimension(:,:,:), allocatable :: x_wrapping_send_buffer, y_wrapping_send_buffer, &
31  x_wrapping_recv_buffer, y_wrapping_recv_buffer

◆ x_wrapping_target_id

integer lowerbc_mod::x_wrapping_target_id
private

Definition at line 34 of file lowerbc.F90.

◆ xx0con

real(kind=default_precision) lowerbc_mod::xx0con
private

Definition at line 27 of file lowerbc.F90.

◆ y2con

real(kind=default_precision) lowerbc_mod::y2con
private

Definition at line 27 of file lowerbc.F90.

◆ y_wrapping_recv_buffer

real(kind=default_precision), dimension(:,:,:), allocatable lowerbc_mod::y_wrapping_recv_buffer
private

Definition at line 30 of file lowerbc.F90.

◆ y_wrapping_send_buffer

real(kind=default_precision), dimension(:,:,:), allocatable lowerbc_mod::y_wrapping_send_buffer
private

Definition at line 30 of file lowerbc.F90.

◆ y_wrapping_target_id

integer lowerbc_mod::y_wrapping_target_id
private

Definition at line 34 of file lowerbc.F90.

◆ yy0con

real(kind=default_precision) lowerbc_mod::yy0con
private

Definition at line 27 of file lowerbc.F90.