MONC
lowerbc.F90
Go to the documentation of this file.
1 
2 module lowerbc_mod
15  use mpi, only: mpi_request_null, mpi_statuses_ignore
16  implicit none
17 
18 #ifndef TEST_MODE
19  private
20 #endif
21 
23 
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
26 
29 
30  real(kind=default_precision), dimension(:,:,:), allocatable :: x_wrapping_send_buffer, y_wrapping_send_buffer, &
32 
33  integer :: iqv ! index for vapour
35 
37 contains
38 
42  lowerbc_get_descriptor%name="lower_bc"
43  lowerbc_get_descriptor%version=0.1
47  end function lowerbc_get_descriptor
48 
49  subroutine initialisation_callback(current_state)
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 
145  end subroutine initialisation_callback
146 
147  subroutine finalisation_callback(current_state)
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)
154  end subroutine finalisation_callback
155 
156  subroutine allocate_applicable_fields(current_state)
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
171  end subroutine allocate_applicable_fields
172 
173  subroutine timestep_callback(current_state)
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
193  end subroutine timestep_callback
194 
195  subroutine compute_lower_boundary_conditions(current_state, current_y_index, current_x_index, zu, zv, zth, th, zq, q)
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
259  end subroutine compute_lower_boundary_conditions
260 
263  subroutine register_async_wrapping_recv_requests(current_state)
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
277 
282  subroutine complete_async_wrapping(current_state, zth, zq)
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
378  end subroutine complete_async_wrapping
379 
386  subroutine package_y_wrapping_send_buffer(current_state, zth, zq, first_y_index, second_y_index)
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
406  end subroutine package_y_wrapping_send_buffer
407 
414  subroutine package_x_wrapping_send_buffer(current_state, zth, zq, first_x_index, second_x_index)
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
434  end subroutine package_x_wrapping_send_buffer
435 
442  subroutine unpackage_y_wrapping_recv_buffer(current_state, zth, zq, first_y_index, second_y_index)
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
462  end subroutine unpackage_y_wrapping_recv_buffer
463 
470  subroutine unpackage_x_wrapping_recv_buffer(current_state, zth, zq, first_x_index, second_x_index)
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
490  end subroutine unpackage_x_wrapping_recv_buffer
491 
492  subroutine handle_convective_fluxes(current_state, current_y_index, current_x_index, horizontal_velocity_at_k2, th, q)
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)= &
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))+th%data(2, current_y_index, current_x_index)-&
514  current_state%global_grid%configuration%vertical%thref(1)+&
515  current_state%global_grid%configuration%vertical%thref(2)
516 
517  ! Surface Flux of vapour
518  if (current_state%number_q_fields .gt. 0) then
519  n=iqv
520  q(n)%data(1, current_y_index, current_x_index)=q(n)%data(2, current_y_index, current_x_index)+&
521  current_state%surface_vapour_flux*current_state%global_grid%configuration%vertical%dzn(2)/&
522  current_state%diff_coefficient%data(1, current_y_index, current_x_index)
523  endif
524  end subroutine handle_convective_fluxes
525 
526  real(kind=default_precision) function look(current_state, vel)
527  type(model_state_type), target, intent(inout) :: current_state
528  real(kind=default_precision), intent(in) :: vel ! Horizontal speed at lowest model level
529 
530  real(kind=default_precision) :: lookup_real_posn
531  integer :: lookup_int_posn
532 
533  lookup_real_posn=1.0_default_precision+real(current_state%lookup_table_entries-1)*&
534  log(vel/current_state%velmin)*current_state%aloginv
535  lookup_int_posn=int(lookup_real_posn)
536 
537  if (lookup_int_posn .ge. 1) then
538  if (lookup_int_posn .lt. current_state%lookup_table_entries) then ! Linear interpolation
539  look=current_state%lookup_table_ustr(lookup_int_posn)+ (lookup_real_posn-real(lookup_int_posn))*&
540  (current_state%lookup_table_ustr(lookup_int_posn+1)-current_state%lookup_table_ustr(lookup_int_posn))
541  else ! Near neutral
542  look=vel*current_state%cneut
543  end if
544  else ! Nearly free convection
545  look=vel**(-convective_limit)*current_state%cfc
546  end if
547  end function look
548 
549 ! check allocations and initialisations
550 
551  subroutine handle_neutral_fluxes(current_state, current_y_index, current_x_index, horizontal_velocity_at_k2, th, q)
552  type(model_state_type), target, intent(inout) :: current_state
553  type(prognostic_field_type), intent(inout) :: th
554  type(prognostic_field_type), dimension(:), intent(inout) :: q
555  integer, intent(in) :: current_y_index, current_x_index
556  real(kind=default_precision) :: horizontal_velocity_at_k2
557 
558  integer :: n
559  real(kind=default_precision) :: ustr
560 
561  ustr=horizontal_velocity_at_k2*current_state%global_grid%configuration%vertical%vk_on_zlogm
562  current_state%vis_coefficient%data(1, current_y_index, current_x_index)=current_state%global_grid%configuration%vertical%czn*&
563  ustr**2/horizontal_velocity_at_k2
564  current_state%diff_coefficient%data(1, current_y_index, current_x_index)=&
565  current_state%vis_coefficient%data(1, current_y_index, current_x_index)*&
566  current_state%global_grid%configuration%vertical%zlogm/(alphah*current_state%global_grid%configuration%vertical%zlogth)
567  if (current_state%th%active) th%data(1, current_y_index, current_x_index)= &
568  (current_state%surface_temperature_flux*current_state%global_grid%configuration%vertical%dzn(2)/&
569  current_state%diff_coefficient%data(1, current_y_index, current_x_index))+th%data(2, current_y_index, current_x_index)-&
570  current_state%global_grid%configuration%vertical%thref(1)+&
571  current_state%global_grid%configuration%vertical%thref(2)
572 
573  ! Flux of vapour
574  if (current_state%number_q_fields .gt. 0) then
575  n=iqv
576  q(n)%data(1, current_y_index, current_x_index)=q(n)%data(2, current_y_index, current_x_index)+&
577  current_state%surface_vapour_flux*current_state%global_grid%configuration%vertical%dzn(2)/&
578  current_state%diff_coefficient%data(1, current_y_index, current_x_index)
579  endif
580  end subroutine handle_neutral_fluxes
581 
582  subroutine handle_stable_fluxes(current_state, current_y_index, current_x_index, horizontal_velocity_at_k2, th, q)
583  type(model_state_type), target, intent(inout) :: current_state
584  type(prognostic_field_type), intent(inout) :: th
585  type(prognostic_field_type), dimension(:), intent(inout) :: q
586  integer, intent(in) :: current_y_index, current_x_index
587  real(kind=default_precision) :: horizontal_velocity_at_k2
588 
589  real(kind=default_precision) :: ustr
590  integer :: n
591 
592  ustr=horizontal_velocity_at_k2*current_state%global_grid%configuration%vertical%vk_on_zlogm
593  if((current_state%fbuoy - 1.e-9_default_precision) .lt. -4.0_default_precision*&
594  von_karman_constant**2*horizontal_velocity_at_k2**3/ (27.0_default_precision*betam*&
595  current_state%global_grid%configuration%vertical%zn(2)*current_state%global_grid%configuration%vertical%zlogm**2)) then
596  ! Too stable for turbulence
597  current_state%vis_coefficient%data(1, current_y_index, current_x_index)=0.0_default_precision
598  current_state%diff_coefficient%data(1, current_y_index, current_x_index)=0.0_default_precision
599  ! Level 1 values of l_th and l_q set to be harmless for advection scheme
600  if (current_state%th%active) th%data(1, current_y_index, current_x_index)=th%data(2, current_y_index, current_x_index)
601  if (current_state%number_q_fields .gt. 0) then
602  do n=1, current_state%number_q_fields
603  q(n)%data(1, current_y_index, current_x_index)=q(n)%data(2, current_y_index, current_x_index)
604  end do
605  end if
606  else
607  ustr=ustr/3.0_default_precision*(1.0_default_precision-2.0_default_precision*&
608  cos((acos(-27.0_default_precision*betam*von_karman_constant*current_state%global_grid%configuration%vertical%zn(2)*&
609  current_state%fbuoy/(current_state%global_grid%configuration%vertical%zlogm*&
610  2.0_default_precision*ustr**3)-1.0_default_precision)+ 2.0_default_precision*pi)/3.0_default_precision))
611  current_state%vis_coefficient%data(1, current_y_index, current_x_index)=&
612  current_state%global_grid%configuration%vertical%czn*ustr**2/horizontal_velocity_at_k2
613  current_state%diff_coefficient%data(1, current_y_index, current_x_index)=&
614  current_state%vis_coefficient%data(1, current_y_index, current_x_index)*&
615  (current_state%global_grid%configuration%vertical%zlogm-betam*current_state%global_grid%configuration%vertical%zn(2)*&
616  von_karman_constant*current_state%fbuoy/ustr**3)/(alphah*current_state%global_grid%configuration%vertical%zlogth-betah*&
617  von_karman_constant*current_state%fbuoy* (current_state%global_grid%configuration%vertical%zn(2)+ z0-z0th)/ustr**3)
618  if (current_state%th%active) th%data(1, current_y_index, current_x_index)= &
619  (current_state%surface_temperature_flux*current_state%global_grid%configuration%vertical%dzn(2)/&
620  current_state%diff_coefficient%data(1, current_y_index, current_x_index))+th%data(2, current_y_index, current_x_index)-&
621  current_state%global_grid%configuration%vertical%thref(1)+&
622  current_state%global_grid%configuration%vertical%thref(2)
623 
624 
625  !Flux of vapour
626  if (current_state%number_q_fields .gt. 0) then
627  n=iqv
628  q(n)%data(1, current_y_index, current_x_index)=q(n)%data(2, current_y_index, current_x_index)+&
629  current_state%surface_vapour_flux*current_state%global_grid%configuration%vertical%dzn(2)/&
630  current_state%diff_coefficient%data(1, current_y_index, current_x_index)
631  endif
632  end if
633  end subroutine handle_stable_fluxes
634 
635  ! set surface_boundary_flux in init == FBUOY
636  subroutine compute_using_fixed_surface_fluxes(current_state, current_y_index, current_x_index, horizontal_velocity_at_k2, th, q)
637  type(model_state_type), target, intent(inout) :: current_state
638  type(prognostic_field_type), intent(inout) :: th
639  type(prognostic_field_type), dimension(:), intent(inout) :: q
640  integer, intent(in) :: current_y_index, current_x_index
641  real(kind=default_precision) :: horizontal_velocity_at_k2
642 
643  if (current_state%fbuoy .gt. 0.0_default_precision) then
644  call handle_convective_fluxes(current_state, current_y_index, current_x_index, horizontal_velocity_at_k2, th, q)
645  else if (current_state%fbuoy .eq. 0.0_default_precision) then
646  call handle_neutral_fluxes(current_state, current_y_index, current_x_index, horizontal_velocity_at_k2, th, q)
647  else
648  call handle_stable_fluxes(current_state, current_y_index, current_x_index, horizontal_velocity_at_k2, th, q)
649  end if
651 
652 
653  subroutine compute_using_fixed_surface_temperature(current_state, current_y_index, current_x_index, horizontal_velocity_at_k2, &
654  zth, th, zq, q)
655  type(model_state_type), target, intent(inout) :: current_state
656  type(prognostic_field_type), intent(inout) :: th, zth
657  type(prognostic_field_type), dimension(:), intent(inout) :: q, zq
658  real(kind=default_precision), intent(in) :: horizontal_velocity_at_k2
659  integer, intent(in) :: current_y_index, current_x_index
660 
661  real(kind=default_precision) :: dthv_surf, ustr, thvstr
662  integer :: convergence_status, n
663 
664  if (current_state%passive_q) then ! i.e. q is not active
665  ! Assuming no liquid water at level 2
666  dthv_surf = zth%data(2, current_y_index, current_x_index) + &
667  current_state%global_grid%configuration%vertical%thref(2) - current_state%theta_virtual_surf
668  else
669  dthv_surf=zth%data(2, current_y_index, current_x_index) + current_state%global_grid%configuration%vertical%thref(2)&
670  *(1.0_default_precision + current_state%cq(current_state%water_vapour_mixing_ratio_index)*&
671  zq(current_state%water_vapour_mixing_ratio_index)%data(2,current_y_index,current_x_index)) - &
672  current_state%theta_virtual_surf
673  end if
674  convergence_status = mostbc(current_state, horizontal_velocity_at_k2, dthv_surf,&
675  current_state%global_grid%configuration%vertical%zn(2), ustr, thvstr)
676 
677  current_state%vis_coefficient%data(1, current_y_index, current_x_index)=&
678  current_state%global_grid%configuration%vertical%czn*ustr**2/horizontal_velocity_at_k2
679  current_state%diff_coefficient%data(1, current_y_index, current_x_index)=&
680  current_state%global_grid%configuration%vertical%czn*ustr*thvstr/(dthv_surf+smallp)
681  zth%data(1, current_y_index, current_x_index)=2.0*current_state%theta_surf-zth%data(2, current_y_index, current_x_index)-&
682  (current_state%global_grid%configuration%vertical%thref(2)+current_state%global_grid%configuration%vertical%thref(1))
683  th%data(1, current_y_index, current_x_index)=2.0*current_state%theta_surf-th%data(2, current_y_index, current_x_index)-&
684  (current_state%global_grid%configuration%vertical%thref(2)+current_state%global_grid%configuration%vertical%thref(1))
685 
686  if (current_state%number_q_fields .gt. 0) then
687  n=iqv
688  zq(n)%data(1, current_y_index, current_x_index)=zq(n)%data(2, current_y_index, current_x_index)
689  q(n)%data(1, current_y_index, current_x_index)=q(n)%data(2, current_y_index, current_x_index)
690  if (.not. current_state%passive_q) then
691  zq(current_state%water_vapour_mixing_ratio_index)%data(1,current_y_index,current_x_index)=&
692  2.0_default_precision*current_state%surface_vapour_mixing_ratio-&
693  zq(current_state%water_vapour_mixing_ratio_index)%data(2,current_y_index,current_x_index)
694  q(current_state%water_vapour_mixing_ratio_index)%data(1,current_y_index,current_x_index)=&
695  2.0_default_precision*current_state%surface_vapour_mixing_ratio-&
696  q(current_state%water_vapour_mixing_ratio_index)%data(2,current_y_index,current_x_index)
697  endif
698  end if
700 
701  subroutine simple_boundary_values(current_state, current_y_index, current_x_index, th, q)
702  type(model_state_type), target, intent(inout) :: current_state
703  type(prognostic_field_type), intent(inout) :: th
704  type(prognostic_field_type), dimension(:), intent(inout) :: q
705  integer, intent(in) :: current_y_index, current_x_index
706 
707  integer :: n
708 
709  current_state%vis_coefficient%data(1, current_y_index, current_x_index)=0.0_default_precision
710  current_state%diff_coefficient%data(1, current_y_index, current_x_index)=0.0_default_precision
711  if(current_state%backscatter) current_state%dis%data(1, current_y_index, current_x_index)=0.0_default_precision
712  if(current_state%backscatter .and. current_state%th%active) &
713  current_state%dis_th%data(1, current_y_index, current_x_index)=0.0_default_precision
714  if (current_state%th%active) then
715  th%data(1, current_y_index, current_x_index) = th%data(2, current_y_index, current_x_index)
716  end if
717  if (current_state%number_q_fields .gt. 0) then
718  do n=1, current_state%number_q_fields
719  q(n)%data(1, current_y_index, current_x_index) = q(n)%data(2, current_y_index, current_x_index)
720  if (current_state%backscatter) current_state%disq(n)%data(1, current_y_index, current_x_index) = 0.0_default_precision
721  end do
722  end if
723  end subroutine simple_boundary_values
724 
740  integer function mostbc(current_state, delu, delt, z1, ustrdg, tstrdg)
741  type(model_state_type), target, intent(inout) :: current_state
742  real(kind=default_precision), intent(in) :: delu, delt, z1
743  real(kind=default_precision), intent(out) :: ustrdg, tstrdg
744 
745  if (delt .lt. 0.0_default_precision) then
746  if (delu .le. smallp) then
747  ustrdg=0.0_default_precision
748  tstrdg=tstrcona*delt
749  else
750  ! The unstable case
751  mostbc=solve_monin_obukhov_unstable_case(delu, delt, current_state%ellmocon, &
752  ustrdg, tstrdg, current_state%global_grid%configuration%vertical)
753  end if
754  else if (delt .gt. 0.0_default_precision) then
755  ! The stable case
756  mostbc=solve_monin_obukhov_stable_case(delu, delt, current_state%global_grid%configuration%vertical%zlogm, &
757  current_state%cmbc, ustrdg, tstrdg)
758  else
759  ! Trivial neutral case
760  ustrdg=current_state%global_grid%configuration%vertical%vk_on_zlogm*delu
761  tstrdg=0.0_default_precision
763  end if
764 
765  if (mostbc .ne. convergence_success) then
767  !call log_log(LOG_WARN, "Richardson number greater than critical value")
768  else if(mostbc .eq. convergence_failure) then
769  call log_log(log_error, "Convergence failure after 200 iterations")
770  end if
771  end if
772  end function mostbc
773 
774  integer function solve_monin_obukhov_unstable_case(delu, delt, ellmocon, ustrdg, tstrdg, vertical_grid)
775  real(kind=default_precision), intent(in) :: delu, delt, ellmocon
776  real(kind=default_precision), intent(out) :: ustrdg, tstrdg
777  type(vertical_grid_configuration_type), intent(inout) :: vertical_grid
778 
779  integer :: i
780  real(kind=default_precision) :: ellmo, psim, psih, x4, xx, xx0, y2, yy, yy0, err_ustr, err_tstr, &
781  ustrl, tstrl, & ! U and T star at start of iteration
782  ustrpl, tstrpl ! U and T star at end of iteration
783 
784  ! First set initial values
785  ustrl=vertical_grid%vk_on_zlogm*delu
786  tstrl=tstrcona*delt
787 
788  ! Now start iteration
789  do i=1, 200
790  ellmo=ustrl*ustrl*ellmocon/tstrl
791 
792  ! Test for possible square root of negative quantity
793  x4=1.0_default_precision-x4con/ellmo
794  if (x4 .lt. 0.0_default_precision) call log_log(log_error, "Negative square root in x4")
795 
796  xx=sqrt(sqrt(x4))
797  xx0=sqrt(sqrt(1.0_default_precision-xx0con / ellmo))
798  psim=2.*( log((xx+1.0_default_precision)/(xx0+1.0_default_precision))-atan(xx)+atan(xx0) )+&
799  log((xx*xx+1.0_default_precision)/(xx0*xx0+1.0_default_precision))
800  ustrpl=von_karman_constant*delu/(vertical_grid%zlogm-psim)
801 
802  ! Test for possible square root of negative quantity
803  y2=1.-y2con/ellmo
804  if (y2 .lt. 0.0_default_precision) call log_log(log_error, "Negative square root in y2")
805  yy=sqrt(y2)
806  yy0=sqrt(1.0_default_precision-yy0con/ellmo)
807  psih=2.*log((1.0_default_precision+yy)/(1.0_default_precision+yy0))
808  tstrpl=tstrconb*delt/(vertical_grid%zlogth-psih)
809  err_ustr=abs((ustrpl-ustrl)/ ustrl)
810  err_tstr=abs((tstrpl-tstrl)/ tstrl)
811  if ((err_tstr .lt. tolt) .and. (err_ustr .lt. tolm)) then
812  ustrdg=ustrpl
813  tstrdg=tstrpl
815  return
816  else
817  ustrl=(1.0_default_precision-smth)*ustrpl+smth*ustrl
818  tstrl=(1.0_default_precision-smth)*tstrpl+smth*tstrl
819  end if
820  end do
823 
824  integer function solve_monin_obukhov_stable_case(delu, delt, zlogm, cmbc, ustrdg, tstrdg)
825  real(kind=default_precision), intent(in) :: delu, delt, zlogm, cmbc
826  real(kind=default_precision), intent(out) :: ustrdg, tstrdg
827 
828  real(kind=default_precision) :: am, ah, ee, ff, det
829 
830  am=von_karman_constant*delu
831  ah=von_karman_constant*delt
832  ee=am*eecon
833  ff=ah*cmbc-rhmbc*am*am !
834  det=ee*ee-ddbc_x4*ff
836  ! Test for laminar flow
837  if (ff .gt. 0.0_default_precision) then
838  if ((ee .lt. 0.0_default_precision).and.(det .gt. 0.0_default_precision)) then
839  ustrdg=(-ee+sqrt(det))*r2ddbc
840  tstrdg=ustrdg*(am-zlogm*ustrdg)*rcmbc
841  else
843  ustrdg=0.0_default_precision
844  tstrdg=0.0_default_precision
845  end if
846  else if (ddbc .eq. 0.0_default_precision) then
847  ! Degenerate case
848  ustrdg=-ff/ee
849  tstrdg=delt*ustrdg/delu
850  else
851  ! Solve quadratic for USTRDG
852  ustrdg=(-ee+sqrt(det))*r2ddbc
853  tstrdg=ustrdg*(am-zlogm*ustrdg)*rcmbc
854  end if
856 end module lowerbc_mod
registry_mod::is_component_enabled
logical function, public is_component_enabled(options_database, component_name)
Determines whether or not a specific component is registered and enabled.
Definition: registry.F90:334
logging_mod::log_error
integer, parameter, public log_error
Only log ERROR messages.
Definition: logging.F90:11
lowerbc_mod::compute_using_fixed_surface_fluxes
subroutine compute_using_fixed_surface_fluxes(current_state, current_y_index, current_x_index, horizontal_velocity_at_k2, th, q)
Definition: lowerbc.F90:637
prognostics_mod
Contains prognostic field definitions and functions.
Definition: prognostics.F90:2
lowerbc_mod::complete_async_wrapping
subroutine complete_async_wrapping(current_state, zth, zq)
Completes the asynchronous wrapping if required for periodic boundary conditions.
Definition: lowerbc.F90:283
logging_mod::log_warn
integer, parameter, public log_warn
Log WARNING and ERROR messages.
Definition: logging.F90:12
state_mod::prescribed_surface_values
integer, parameter, public prescribed_surface_values
Definition: state.F90:15
lowerbc_mod
This sets the lower boundary conditions for theta and the q variables.
Definition: lowerbc.F90:2
lowerbc_mod::wrapping_comm_requests
integer, dimension(4) wrapping_comm_requests
Definition: lowerbc.F90:34
lowerbc_mod::tstrcona
real(kind=default_precision) tstrcona
Definition: lowerbc.F90:27
lowerbc_mod::handle_convective_fluxes
subroutine handle_convective_fluxes(current_state, current_y_index, current_x_index, horizontal_velocity_at_k2, th, q)
Definition: lowerbc.F90:493
grids_mod::x_index
integer, parameter, public x_index
Definition: grids.F90:14
grids_mod::y_index
integer, parameter, public y_index
Definition: grids.F90:14
lowerbc_mod::y_wrapping_send_buffer
real(kind=default_precision), dimension(:,:,:), allocatable y_wrapping_send_buffer
Definition: lowerbc.F90:30
science_constants_mod::gammam
real(kind=default_precision), public gammam
Definition: scienceconstants.F90:13
lowerbc_mod::finalisation_callback
subroutine finalisation_callback(current_state)
Definition: lowerbc.F90:148
datadefn_mod::precision_type
integer, public precision_type
Definition: datadefn.F90:19
lowerbc_mod::y2con
real(kind=default_precision) y2con
Definition: lowerbc.F90:27
lowerbc_mod::timestep_callback
subroutine timestep_callback(current_state)
Definition: lowerbc.F90:174
logging_mod::log_log
subroutine, public log_log(level, message, str)
Logs a message at the specified level. If the level is above the current level then the message is ig...
Definition: logging.F90:75
lowerbc_mod::package_x_wrapping_send_buffer
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.
Definition: lowerbc.F90:415
lowerbc_mod::iqv
integer iqv
Definition: lowerbc.F90:33
lowerbc_mod::x_wrapping_send_buffer
real(kind=default_precision), dimension(:,:,:), allocatable x_wrapping_send_buffer
Definition: lowerbc.F90:30
lowerbc_mod::tstrconb
real(kind=default_precision) tstrconb
Definition: lowerbc.F90:27
lowerbc_mod::x4con
real(kind=default_precision) x4con
Definition: lowerbc.F90:27
lowerbc_mod::ddbc_x4
real(kind=default_precision) ddbc_x4
Definition: lowerbc.F90:27
lowerbc_mod::register_async_wrapping_recv_requests
subroutine register_async_wrapping_recv_requests(current_state)
Registers asynchronous wrapping recv requests as needed.
Definition: lowerbc.F90:264
lowerbc_mod::unpackage_x_wrapping_recv_buffer
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.
Definition: lowerbc.F90:471
lowerbc_mod::compute_using_fixed_surface_temperature
subroutine compute_using_fixed_surface_temperature(current_state, current_y_index, current_x_index, horizontal_velocity_at_k2, zth, th, zq, q)
Definition: lowerbc.F90:655
monc_component_mod
Interfaces and types that MONC components must specify.
Definition: monc_component.F90:6
lowerbc_mod::rhmbc
real(kind=default_precision) rhmbc
Definition: lowerbc.F90:27
lowerbc_mod::viscous_courant_coefficient
real(kind=default_precision) viscous_courant_coefficient
Definition: lowerbc.F90:27
lowerbc_mod::unpackage_y_wrapping_recv_buffer
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.
Definition: lowerbc.F90:443
lowerbc_mod::tolt
real(kind=default_precision), parameter tolt
Definition: lowerbc.F90:24
lowerbc_mod::tolm
real(kind=default_precision), parameter tolm
Definition: lowerbc.F90:24
lowerbc_mod::y_wrapping_recv_buffer
real(kind=default_precision), dimension(:,:,:), allocatable y_wrapping_recv_buffer
Definition: lowerbc.F90:30
lowerbc_mod::handle_stable_fluxes
subroutine handle_stable_fluxes(current_state, current_y_index, current_x_index, horizontal_velocity_at_k2, th, q)
Definition: lowerbc.F90:583
lowerbc_mod::solve_monin_obukhov_unstable_case
integer function solve_monin_obukhov_unstable_case(delu, delt, ellmocon, ustrdg, tstrdg, vertical_grid)
Definition: lowerbc.F90:775
science_constants_mod
Scientific constant values used throughout simulations. Each has a default value and this can be over...
Definition: scienceconstants.F90:3
lowerbc_mod::ddbc
real(kind=default_precision) ddbc
Definition: lowerbc.F90:27
grids_mod::z_index
integer, parameter, public z_index
Grid index parameters.
Definition: grids.F90:14
science_constants_mod::smallp
real(kind=default_precision), public smallp
Definition: scienceconstants.F90:13
q_indices_mod::standard_q_names
type(standard_q_names_type), public standard_q_names
Definition: q_indices.F90:59
lowerbc_mod::rcmbc
real(kind=default_precision) rcmbc
Definition: lowerbc.F90:27
state_mod::model_state_type
The ModelState which represents the current state of a run.
Definition: state.F90:39
state_mod::forward_stepping
integer, parameter, public forward_stepping
Definition: state.F90:15
lowerbc_mod::convergence_failure
integer, parameter convergence_failure
Definition: lowerbc.F90:22
lowerbc_mod::xx0con
real(kind=default_precision) xx0con
Definition: lowerbc.F90:27
lowerbc_mod::smth
real(kind=default_precision), parameter smth
Definition: lowerbc.F90:24
lowerbc_mod::allocate_applicable_fields
subroutine allocate_applicable_fields(current_state)
Definition: lowerbc.F90:157
lowerbc_mod::x_wrapping_target_id
integer x_wrapping_target_id
Definition: lowerbc.F90:34
lowerbc_mod::lowerbc_get_descriptor
type(component_descriptor_type) function, public lowerbc_get_descriptor()
Descriptor of this component for registration.
Definition: lowerbc.F90:42
lowerbc_mod::y_wrapping_target_id
integer y_wrapping_target_id
Definition: lowerbc.F90:34
lowerbc_mod::handle_neutral_fluxes
subroutine handle_neutral_fluxes(current_state, current_y_index, current_x_index, horizontal_velocity_at_k2, th, q)
Definition: lowerbc.F90:552
science_constants_mod::betah
real(kind=default_precision), public betah
Definition: scienceconstants.F90:13
q_indices_mod::get_q_index
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
science_constants_mod::betam
real(kind=default_precision), public betam
Definition: scienceconstants.F90:13
science_constants_mod::convective_limit
real(kind=default_precision), public convective_limit
Definition: scienceconstants.F90:13
logging_mod
Logging utility.
Definition: logging.F90:2
prognostics_mod::prognostic_field_type
A prognostic field which is assumed to be 3D.
Definition: prognostics.F90:13
grids_mod::vertical_grid_configuration_type
The configuration of the grid vertically.
Definition: grids.F90:28
datadefn_mod
Contains common definitions for the data and datatypes used by MONC.
Definition: datadefn.F90:2
state_mod::prescribed_surface_fluxes
integer, parameter, public prescribed_surface_fluxes
Definition: state.F90:15
science_constants_mod::z0th
real(kind=default_precision), public z0th
Definition: scienceconstants.F90:13
lowerbc_mod::compute_lower_boundary_conditions
subroutine compute_lower_boundary_conditions(current_state, current_y_index, current_x_index, zu, zv, zth, th, zq, q)
Definition: lowerbc.F90:196
logging_mod::log_master_log
subroutine, public log_master_log(level, message)
Will log just from the master process.
Definition: logging.F90:47
q_indices_mod
This manages the Q variables and specifically the mapping between names and the index that they are s...
Definition: q_indices.F90:2
lowerbc_mod::eecon
real(kind=default_precision) eecon
Definition: lowerbc.F90:27
lowerbc_mod::yy0con
real(kind=default_precision) yy0con
Definition: lowerbc.F90:27
registry_mod
MONC component registry.
Definition: registry.F90:5
science_constants_mod::alphah
real(kind=default_precision), public alphah
Definition: scienceconstants.F90:13
lowerbc_mod::simple_boundary_values
subroutine simple_boundary_values(current_state, current_y_index, current_x_index, th, q)
Definition: lowerbc.F90:702
lowerbc_mod::mostbc
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 ...
Definition: lowerbc.F90:741
lowerbc_mod::convergence_richardson_too_large
integer, parameter convergence_richardson_too_large
Definition: lowerbc.F90:22
science_constants_mod::z0
real(kind=default_precision), public z0
Definition: scienceconstants.F90:13
lowerbc_mod::r2ddbc
real(kind=default_precision) r2ddbc
Definition: lowerbc.F90:27
lowerbc_mod::x_wrapping_recv_buffer
real(kind=default_precision), dimension(:,:,:), allocatable x_wrapping_recv_buffer
Definition: lowerbc.F90:30
lowerbc_mod::convergence_success
integer, parameter convergence_success
Definition: lowerbc.F90:22
grids_mod
Functionality to support the different types of grid and abstraction between global grids and local o...
Definition: grids.F90:5
science_constants_mod::von_karman_constant
real(kind=default_precision), public von_karman_constant
Definition: scienceconstants.F90:13
lowerbc_mod::solve_monin_obukhov_stable_case
integer function solve_monin_obukhov_stable_case(delu, delt, zlogm, cmbc, ustrdg, tstrdg)
Definition: lowerbc.F90:825
lowerbc_mod::initialisation_callback
subroutine initialisation_callback(current_state)
Definition: lowerbc.F90:50
lowerbc_mod::package_y_wrapping_send_buffer
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.
Definition: lowerbc.F90:387
science_constants_mod::pi
real(kind=default_precision), public pi
Definition: scienceconstants.F90:13
science_constants_mod::gammah
real(kind=default_precision), public gammah
Definition: scienceconstants.F90:13
lowerbc_mod::look
real(kind=default_precision) function look(current_state, vel)
Definition: lowerbc.F90:527
monc_component_mod::component_descriptor_type
Description of a component.
Definition: monc_component.F90:42
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
state_mod
The model state which represents the current state of a run.
Definition: state.F90:2