MONC
Functions/Subroutines | Variables
smagorinsky_mod Module Reference

Calculates the Smagorinsky eddy viscosity and diffusivity at l_w-points. More...

Functions/Subroutines

type(component_descriptor_type) function, public smagorinsky_get_descriptor ()
 Provides the descriptor back to the caller and is used in component registration. More...
 
subroutine initialisation_callback (current_state)
 Initialisation call back which will read in the coriolis configuration and set up the geostrophic winds. More...
 
subroutine timestep_callback (current_state)
 For each none halo cell this will calculate the subgrid terms for viscosity and diffusion. More...
 
subroutine finalisation_callback (current_state)
 Called when the model is finishing up, will finalise the halo communications represented by the state. More...
 
subroutine update_viscous_number (current_state)
 Update viscous number based upon the viscosity and diffusivity coefficients. More...
 
subroutine setfri (current_state, richardson_number, ssq)
 Calculates the eddy viscosity (VIS) and diffusivity (DIFF) depending on the Richardson Number (RI) and half squared strain rate. More...
 
real(kind=default_precision) function, dimension(current_state%local_grid%size(z_index)), public calculate_richardson_number (current_state, ssq, th, q)
 Calculates the richardson number depending upon the setup of the model and the method selected. More...
 
real(kind=default_precision) function, dimension(current_state%local_grid%size(z_index)) moist_ri_2 (current_state, ssq, th, q)
 Calculates another "moist version" of the Richardson number based on "change in subgrid buoyancy flux when a fraction EPS is exchanged". More...
 
real(kind=default_precision) function, dimension(current_state%local_grid%size(z_index)) moist_ri_1 (current_state, ssq, th, q)
 Calculates numerator of Richarson number as the difference in buoyancy between a parcel at K+1 and a parcel lifted from K to K+1. Condensation is calculated by assuming that TL=T-(Lv/CP)*QL+gz/CP and QT = QV+QL are conserved on lifting. Microphysics conversions are assumed to happen slowly compared to turbulent processes so are, therefore, neglected in the calculation. More...
 
real(kind=default_precision) function, dimension(current_state%local_grid%size(z_index)), public calculate_half_squared_strain_rate (current_state, u, v, w)
 Calculates the half squared strain rate on l_w-points which is used in determining the Richardson number. CX=1./DX, CY=1./DY, RDZ(K)=1./DZ(K), RDZN(K) =1./DZN(K) _SSQ= 0.5*^^DU_I/DX_J+DU_J/DX_I^^**2 _SSQIJ= (DU_I/DX_J+DU_J/DX_I)**2 _Hence SSQ= SUM(I,J) {0.5*(SSQIJ)}. More...
 
real(kind=default_precision) function, dimension(current_state%local_grid%size(z_index)), public calculate_thermal_dissipation_rate (current_state, th)
 
subroutine copy_vis_to_halo_buffer (current_state, neighbour_description, dim, source_index, pid_location, current_page, source_data)
 Copies the viscosity field data to halo buffers for a specific process in a dimension and halo cell. More...
 
subroutine copy_vis_corners_to_halo_buffer (current_state, neighbour_description, corner_loc, x_source_index, y_source_index, pid_location, current_page, source_data)
 Copies the viscosity field corner data to halo buffers for a specific process. More...
 
subroutine copy_diff_to_halo_buffer (current_state, neighbour_description, dim, source_index, pid_location, current_page, source_data)
 Copies the diffusion field data to halo buffers for a specific process in a dimension and halo cell. More...
 
subroutine copy_diff_corners_to_halo_buffer (current_state, neighbour_description, corner_loc, x_source_index, y_source_index, pid_location, current_page, source_data)
 Copies the diffusion field corner data to halo buffers for a specific process. More...
 

Variables

integer, parameter richardson_number_calculation =2
 
real(kind=default_precision) eps
 
real(kind=default_precision) repsh
 
real(kind=default_precision) thcona
 
real(kind=default_precision) thconb
 
real(kind=default_precision) thconap1
 
real(kind=default_precision) suba
 
real(kind=default_precision) subb
 
real(kind=default_precision) subc
 
real(kind=default_precision) subg
 
real(kind=default_precision) subh
 
real(kind=default_precision) subr
 
real(kind=default_precision) pr_n
 
real(kind=default_precision) ric
 
real(kind=default_precision) ricinv
 
integer iqv
 
integer iql
 

Detailed Description

Calculates the Smagorinsky eddy viscosity and diffusivity at l_w-points.

Function/Subroutine Documentation

◆ calculate_half_squared_strain_rate()

real(kind=default_precision) function, dimension(current_state%local_grid%size(z_index)), public smagorinsky_mod::calculate_half_squared_strain_rate ( type(model_state_type), intent(inout), target  current_state,
type(prognostic_field_type), intent(inout)  u,
type(prognostic_field_type), intent(inout)  v,
type(prognostic_field_type), intent(inout)  w 
)

Calculates the half squared strain rate on l_w-points which is used in determining the Richardson number. CX=1./DX, CY=1./DY, RDZ(K)=1./DZ(K), RDZN(K) =1./DZN(K) _SSQ= 0.5*^^DU_I/DX_J+DU_J/DX_I^^**2 _SSQIJ= (DU_I/DX_J+DU_J/DX_I)**2 _Hence SSQ= SUM(I,J) {0.5*(SSQIJ)}.

Parameters
current_stateThe current model state
Returns
The hard squared strain rate for a specific column

Definition at line 410 of file smagorinsky.F90.

411  type(model_state_type), target, intent(inout) :: current_state
412  type(prognostic_field_type), intent(inout) :: u, v, w
413  real(kind=default_precision) :: calculate_half_squared_strain_rate(current_state%local_grid%size(z_index))
414 
415  integer :: k
416  real(kind=default_precision) :: ssq11, ssq22, ssq33, ssq13, ssq23, ssq12
417 
418  do k=2,current_state%local_grid%size(z_index)-1
419 #ifdef U_ACTIVE
420  ssq11=current_state%global_grid%configuration%horizontal%cx2*(&
421  (u%data(k+1, current_state%column_local_y, current_state%column_local_x)-&
422  u%data(k+1, current_state%column_local_y, current_state%column_local_x-1))**2+&
423  (u%data(k, current_state%column_local_y, current_state%column_local_x)-&
424  u%data(k, current_state%column_local_y, current_state%column_local_x-1))**2)
425 #else
426  ssq11=0.0_default_precision
427 #endif
428 #ifdef V_ACTIVE
429  ssq22=current_state%global_grid%configuration%horizontal%cy2*(&
430  (v%data(k+1, current_state%column_local_y, current_state%column_local_x)-&
431  v%data(k+1, current_state%column_local_y-1, current_state%column_local_x))**2+&
432  (v%data(k, current_state%column_local_y, current_state%column_local_x)-&
433  v%data(k, current_state%column_local_y-1, current_state%column_local_x))**2)
434 #else
435  ssq22=0.0_default_precision
436 #endif
437 #ifdef W_ACTIVE
438  ssq33=((w%data(k, current_state%column_local_y, current_state%column_local_x)-&
439  w%data(k-1, current_state%column_local_y, current_state%column_local_x))*&
440  current_state%global_grid%configuration%vertical%rdz(k))**2 +&
441  ((w%data(k+1, current_state%column_local_y, current_state%column_local_x)-&
442  w%data(k, current_state%column_local_y, current_state%column_local_x))*&
443  current_state%global_grid%configuration%vertical%rdz(k+1))**2
444 #else
445  ssq33=0.0_default_precision
446 #endif
447 #if defined(U_ACTIVE) && defined(W_ACTIVE)
448  ! Average over 2 points U and W
449  ssq13=(((u%data(k+1, current_state%column_local_y, current_state%column_local_x)-&
450  u%data(k, current_state%column_local_y, current_state%column_local_x))*&
451  current_state%global_grid%configuration%vertical%rdzn(k+1)+&
452  (w%data(k, current_state%column_local_y, current_state%column_local_x+1)-&
453  w%data(k, current_state%column_local_y, current_state%column_local_x))*&
454  current_state%global_grid%configuration%horizontal%cx)**2+&
455  ((u%data(k+1, current_state%column_local_y, current_state%column_local_x-1)-&
456  u%data(k, current_state%column_local_y, current_state%column_local_x-1))*&
457  current_state%global_grid%configuration%vertical%rdzn(k+1)+&
458  (w%data(k, current_state%column_local_y, current_state%column_local_x)-&
459  w%data(k, current_state%column_local_y, current_state%column_local_x-1))*&
460  current_state%global_grid%configuration%horizontal%cx)**2)*0.5_default_precision
461 #elif defined(U_ACTIVE) && !defined(W_ACTIVE)
462  ssq13=(((u%data(k+1, current_state%column_local_y, current_state%column_local_x)-&
463  u%data(k, current_state%column_local_y, current_state%column_local_x))*&
464  current_state%global_grid%configuration%vertical%rdzn(k+1))**2+&
465  ((u%data(k+1, current_state%column_local_y, current_state%column_local_x-1)-&
466  u%data(k, current_state%column_local_y, current_state%column_local_x-1))*&
467  current_state%global_grid%configuration%vertical%rdzn(k+1))**2)*0.5_default_precision
468 #elif !defined(U_ACTIVE) && defined(W_ACTIVE)
469  ssq13=(((w%data(k, current_state%column_local_y, current_state%column_local_x+1)-&
470  w%data(k, current_state%column_local_y, current_state%column_local_x))*&
471  current_state%global_grid%configuration%horizontal%cx)**2+&
472  ((w%data(k, current_state%column_local_y, current_state%column_local_x)-&
473  w%data(k, current_state%column_local_y, current_state%column_local_x-1))*&
474  current_state%global_grid%configuration%horizontal%cx)**2)*0.5_default_precision
475 #else
476  ssq13=0.0_default_precision
477 #endif
478 #if defined(W_ACTIVE) && defined(V_ACTIVE)
479  ! Average over 2 points W and V
480  ssq23=(((w%data(k, current_state%column_local_y, current_state%column_local_x)-&
481  w%data(k, current_state%column_local_y-1, current_state%column_local_x))*&
482  current_state%global_grid%configuration%horizontal%cy+&
483  (v%data(k+1, current_state%column_local_y-1, current_state%column_local_x)-&
484  v%data(k, current_state%column_local_y-1, current_state%column_local_x))*&
485  current_state%global_grid%configuration%vertical%rdzn(k+1))**2+&
486  ((w%data(k, current_state%column_local_y+1, current_state%column_local_x)-&
487  w%data(k, current_state%column_local_y, current_state%column_local_x))*&
488  current_state%global_grid%configuration%horizontal%cy+&
489  (v%data(k+1, current_state%column_local_y, current_state%column_local_x)-&
490  v%data(k, current_state%column_local_y, current_state%column_local_x))*&
491  current_state%global_grid%configuration%vertical%rdzn(k+1))**2)*0.5_default_precision
492 #elif defined(W_ACTIVE) && !defined(V_ACTIVE)
493  ssq23=(((w%data(k, current_state%column_local_y, current_state%column_local_x)-&
494  w%data(k, current_state%column_local_y-1, current_state%column_local_x))*&
495  current_state%global_grid%configuration%horizontal%cy)**2+&
496  ((w%data(k, current_state%column_local_y+1, current_state%column_local_x)-&
497  w%data(k, current_state%column_local_y, current_state%column_local_x))*&
498  current_state%global_grid%configuration%horizontal%cy)**2)*0.5_default_precision
499 #elif !defined(W_ACTIVE) && defined(V_ACTIVE)
500  ssq23=(((v%data(k+1, current_state%column_local_y-1, current_state%column_local_x)-&
501  v%data(k, current_state%column_local_y-1, current_state%column_local_x))*&
502  current_state%global_grid%configuration%vertical%rdzn(k+1))**2+&
503  ((v%data(k+1, current_state%column_local_y, current_state%column_local_x)-&
504  v%data(k, current_state%column_local_y, current_state%column_local_x))*&
505  current_state%global_grid%configuration%vertical%rdzn(k+1))**2)*0.5_default_precision
506 #else
507  ssq23=0.0_default_precision
508 #endif
509 
510 #if defined(U_ACTIVE) && defined(V_ACTIVE)
511  ! Average over 8 points from U and V
512  ssq12=(((((u%data(k, current_state%column_local_y, current_state%column_local_x-1)-&
513  u%data(k, current_state%column_local_y-1, current_state%column_local_x-1))*&
514  current_state%global_grid%configuration%horizontal%cy+&
515  (v%data(k, current_state%column_local_y-1, current_state%column_local_x)-&
516  v%data(k, current_state%column_local_y-1, current_state%column_local_x-1))*&
517  current_state%global_grid%configuration%horizontal%cx)**2 +&
518  ((u%data(k, current_state%column_local_y+1, current_state%column_local_x-1)-&
519  u%data(k, current_state%column_local_y, current_state%column_local_x-1))*&
520  current_state%global_grid%configuration%horizontal%cy+&
521  (v%data(k, current_state%column_local_y, current_state%column_local_x)-&
522  v%data(k, current_state%column_local_y, current_state%column_local_x-1))*&
523  current_state%global_grid%configuration%horizontal%cx)**2) +(&
524  ((u%data(k, current_state%column_local_y, current_state%column_local_x)-&
525  u%data(k, current_state%column_local_y-1, current_state%column_local_x))*&
526  current_state%global_grid%configuration%horizontal%cy+&
527  (v%data(k, current_state%column_local_y-1, current_state%column_local_x+1)-&
528  v%data(k, current_state%column_local_y-1, current_state%column_local_x))*&
529  current_state%global_grid%configuration%horizontal%cx)**2 +&
530  ((u%data(k, current_state%column_local_y+1, current_state%column_local_x)-&
531  u%data(k, current_state%column_local_y, current_state%column_local_x))*&
532  current_state%global_grid%configuration%horizontal%cy+&
533  (v%data(k, current_state%column_local_y, current_state%column_local_x+1)-&
534  v%data(k, current_state%column_local_y, current_state%column_local_x))*&
535  current_state%global_grid%configuration%horizontal%cx)**2))+((&
536  ((u%data(k+1, current_state%column_local_y, current_state%column_local_x-1)-&
537  u%data(k+1, current_state%column_local_y-1, current_state%column_local_x-1))*&
538  current_state%global_grid%configuration%horizontal%cy+&
539  (v%data(k+1, current_state%column_local_y-1, current_state%column_local_x)-&
540  v%data(k+1, current_state%column_local_y-1, current_state%column_local_x-1))*&
541  current_state%global_grid%configuration%horizontal%cx)**2+&
542  ((u%data(k+1, current_state%column_local_y+1, current_state%column_local_x-1)-&
543  u%data(k+1, current_state%column_local_y, current_state%column_local_x-1))*&
544  current_state%global_grid%configuration%horizontal%cy+&
545  (v%data(k+1, current_state%column_local_y, current_state%column_local_x)-&
546  v%data(k+1, current_state%column_local_y, current_state%column_local_x-1))*&
547  current_state%global_grid%configuration%horizontal%cx)**2)+(&
548  ((u%data(k+1, current_state%column_local_y, current_state%column_local_x)-&
549  u%data(k+1, current_state%column_local_y-1, current_state%column_local_x))*&
550  current_state%global_grid%configuration%horizontal%cy+&
551  (v%data(k+1, current_state%column_local_y-1, current_state%column_local_x+1)-&
552  v%data(k+1, current_state%column_local_y-1, current_state%column_local_x))*&
553  current_state%global_grid%configuration%horizontal%cx)**2+&
554  ((u%data(k+1, current_state%column_local_y+1, current_state%column_local_x)-&
555  u%data(k+1, current_state%column_local_y, current_state%column_local_x))*&
556  current_state%global_grid%configuration%horizontal%cy+&
557  (v%data(k+1, current_state%column_local_y, current_state%column_local_x+1)-&
558  v%data(k+1, current_state%column_local_y, current_state%column_local_x))*&
559  current_state%global_grid%configuration%horizontal%cx)**2)))*0.125_default_precision
560 
561 #elif defined(U_ACTIVE) && !defined(V_ACTIVE)
562 
563  ssq12=(((((u%data(k, current_state%column_local_y, current_state%column_local_x-1)-&
564  u%data(k, current_state%column_local_y-1, current_state%column_local_x-1))*&
565  current_state%global_grid%configuration%horizontal%cy)**2 +&
566  ((u%data(k, current_state%column_local_y+1, current_state%column_local_x-1)-&
567  u%data(k, current_state%column_local_y, current_state%column_local_x-1))*&
568  current_state%global_grid%configuration%horizontal%cy)**2) +(&
569  ((u%data(k, current_state%column_local_y, current_state%column_local_x)-&
570  u%data(k, current_state%column_local_y-1, current_state%column_local_x))*&
571  current_state%global_grid%configuration%horizontal%cy)**2 +&
572  ((u%data(k, current_state%column_local_y+1, current_state%column_local_x)-&
573  u%data(k, current_state%column_local_y, current_state%column_local_x))*&
574  current_state%global_grid%configuration%horizontal%cy)**2))+((&
575  ((u%data(k+1, current_state%column_local_y, current_state%column_local_x-1)-&
576  u%data(k+1, current_state%column_local_y-1, current_state%column_local_x-1))*&
577  current_state%global_grid%configuration%horizontal%cy)**2+&
578  ((u%data(k+1, current_state%column_local_y+1, current_state%column_local_x-1)-&
579  u%data(k+1, current_state%column_local_y, current_state%column_local_x-1))*&
580  current_state%global_grid%configuration%horizontal%cy)**2)+(&
581  ((u%data(k+1, current_state%column_local_y, current_state%column_local_x)-&
582  u%data(k+1, current_state%column_local_y-1, current_state%column_local_x))*&
583  current_state%global_grid%configuration%horizontal%cy)**2+&
584  ((u%data(k+1, current_state%column_local_y+1, current_state%column_local_x)-&
585  u%data(k+1, current_state%column_local_y, current_state%column_local_x))*&
586  current_state%global_grid%configuration%horizontal%cy)**2)))*0.125_default_precision
587 
588 #elif !defined(U_ACTIVE) && defined(V_ACTIVE)
589 
590  ssq12=(((((v%data(k, current_state%column_local_y-1, current_state%column_local_x)-&
591  v%data(k, current_state%column_local_y-1, current_state%column_local_x-1))*&
592  current_state%global_grid%configuration%horizontal%cx)**2 +&
593  ((v%data(k, current_state%column_local_y, current_state%column_local_x)-&
594  v%data(k, current_state%column_local_y, current_state%column_local_x-1))*&
595  current_state%global_grid%configuration%horizontal%cx)**2) +&
596  (((v%data(k, current_state%column_local_y-1, current_state%column_local_x+1)-&
597  v%data(k, current_state%column_local_y-1, current_state%column_local_x))*&
598  current_state%global_grid%configuration%horizontal%cx)**2 +&
599  ((v%data(k, current_state%column_local_y, current_state%column_local_x+1)-&
600  v%data(k, current_state%column_local_y, current_state%column_local_x))*&
601  current_state%global_grid%configuration%horizontal%cx)**2))+((&
602  ((v%data(k+1, current_state%column_local_y-1, current_state%column_local_x)-&
603  v%data(k+1, current_state%column_local_y-1, current_state%column_local_x-1))*&
604  current_state%global_grid%configuration%horizontal%cx)**2+&
605  ((v%data(k+1, current_state%column_local_y, current_state%column_local_x)-&
606  v%data(k+1, current_state%column_local_y, current_state%column_local_x-1))*&
607  current_state%global_grid%configuration%horizontal%cx)**2)+(&
608  ((v%data(k+1, current_state%column_local_y-1, current_state%column_local_x+1)-&
609  v%data(k+1, current_state%column_local_y-1, current_state%column_local_x))*&
610  current_state%global_grid%configuration%horizontal%cx)**2+&
611  ((v%data(k+1, current_state%column_local_y, current_state%column_local_x+1)-&
612  v%data(k+1, current_state%column_local_y, current_state%column_local_x))*&
613  current_state%global_grid%configuration%horizontal%cx)**2)))*0.125_default_precision
614 #else
615  ssq12=0.0_default_precision
616 #endif
617  calculate_half_squared_strain_rate(k)=ssq11+ssq22+ssq33+ssq13+ssq23+ssq12+smallp
618  end do
Here is the caller graph for this function:

◆ calculate_richardson_number()

real(kind=default_precision) function, dimension(current_state%local_grid%size(z_index)), public smagorinsky_mod::calculate_richardson_number ( type(model_state_type), intent(inout), target  current_state,
real(kind=default_precision), dimension(:), intent(in)  ssq,
type(prognostic_field_type), intent(inout)  th,
type(prognostic_field_type), dimension(:), intent(inout)  q 
)

Calculates the richardson number depending upon the setup of the model and the method selected.

Parameters
current_stateThe current model state
ssqThe half squared strain rate
Returns
The Richardson number for a specific column

Definition at line 204 of file smagorinsky.F90.

205  type(model_state_type), target, intent(inout) :: current_state
206  real(kind=default_precision), dimension(:), intent(in) :: ssq
207  type(prognostic_field_type), intent(inout) :: th
208  type(prognostic_field_type), dimension(:), intent(inout) :: q
209  real(kind=default_precision), dimension(current_state%local_grid%size(Z_INDEX)) :: calculate_richardson_number
210 
211  integer :: k
212 
213  if (.not. current_state%passive_th) then
214  if (.not. current_state%passive_q) then
215  if (richardson_number_calculation .eq. 2) then
216  calculate_richardson_number=moist_ri_2(current_state, ssq, th, q)
217  else
218  calculate_richardson_number=moist_ri_1(current_state, ssq, th, q)
219  end if
220  else
221  do k=2, current_state%local_grid%size(z_index)-1
222  calculate_richardson_number(k) = (current_state%global_grid%configuration%vertical%dthref(k) + &
223  current_state%th%data(k+1, current_state%column_local_y, current_state%column_local_x) -&
224  current_state%th%data(k, current_state%column_local_y, current_state%column_local_x))* &
225  current_state%global_grid%configuration%vertical%rdzn(k+1)*&
226  current_state%global_grid%configuration%vertical%buoy_co(k) / ssq(k)
227  end do
228  end if
229  else
230  calculate_richardson_number=0.0_default_precision
231  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ calculate_thermal_dissipation_rate()

real(kind=default_precision) function, dimension(current_state%local_grid%size(z_index)), public smagorinsky_mod::calculate_thermal_dissipation_rate ( type(model_state_type), intent(inout), target  current_state,
type(prognostic_field_type), intent(inout)  th 
)

Definition at line 624 of file smagorinsky.F90.

625  type(model_state_type), target, intent(inout) :: current_state
626  type(prognostic_field_type), intent(inout) :: th
627  real(kind=default_precision) :: calculate_thermal_dissipation_rate(current_state%local_grid%size(z_index))
628 
629  integer :: k
630 
631  do k=2,current_state%local_grid%size(z_index)-1
632  calculate_thermal_dissipation_rate(k) = &
633  current_state%diff_coefficient%data(k, current_state%column_local_y, current_state%column_local_x) * ( &
634  (current_state%global_grid%configuration%vertical%rdzn(k+1) * &
635  (current_state%th%data(k+1, current_state%column_local_y, current_state%column_local_x) - &
636  current_state%th%data(k, current_state%column_local_y, current_state%column_local_x) - &
637  current_state%global_grid%configuration%vertical%dthref(k) ) ) ** 2 + &
638  0.25 * current_state%global_grid%configuration%horizontal%cx2 * &
639  ( (current_state%th%data(k, current_state%column_local_y, current_state%column_local_x+1) - &
640  current_state%th%data(k, current_state%column_local_y, current_state%column_local_x) ) ** 2 + &
641  (current_state%th%data(k, current_state%column_local_y, current_state%column_local_x) - &
642  current_state%th%data(k, current_state%column_local_y, current_state%column_local_x-1) ) ** 2 + &
643  (current_state%th%data(k+1, current_state%column_local_y, current_state%column_local_x+1) - &
644  current_state%th%data(k+1, current_state%column_local_y, current_state%column_local_x) ) **2 + &
645  (current_state%th%data(k+1, current_state%column_local_y, current_state%column_local_x) - &
646  current_state%th%data(k+1, current_state%column_local_y, current_state%column_local_x-1) ) ** 2 ) + &
647  0.25 * current_state%global_grid%configuration%horizontal%cy2 * &
648  ( (current_state%th%data(k, current_state%column_local_y+1, current_state%column_local_x) - &
649  current_state%th%data(k, current_state%column_local_y, current_state%column_local_x) ) ** 2 + &
650  (current_state%th%data(k, current_state%column_local_y, current_state%column_local_x) - &
651  current_state%th%data(k, current_state%column_local_y-1, current_state%column_local_x) ) ** 2 + &
652  (current_state%th%data(k+1, current_state%column_local_y+1, current_state%column_local_x) - &
653  current_state%th%data(k+1, current_state%column_local_y, current_state%column_local_x) ) **2 + &
654  (current_state%th%data(k+1, current_state%column_local_y, current_state%column_local_x) - &
655  current_state%th%data(k+1, current_state%column_local_y-1, current_state%column_local_x) ) ** 2 ) )
656  end do
657 

◆ copy_diff_corners_to_halo_buffer()

subroutine smagorinsky_mod::copy_diff_corners_to_halo_buffer ( type(model_state_type), intent(inout)  current_state,
type(neighbour_description_type), intent(inout)  neighbour_description,
integer, intent(in)  corner_loc,
integer, intent(in)  x_source_index,
integer, intent(in)  y_source_index,
integer, intent(in)  pid_location,
integer, dimension(:), intent(inout)  current_page,
type(field_data_wrapper_type), dimension(:), intent(in), optional  source_data 
)
private

Copies the diffusion field corner data to halo buffers for a specific process.

Parameters
current_stateThe current model state
neighbour_descriptionsDescription of the neighbour halo swapping status
corner_locThe corner location
x_source_indexThe X source index of the dimension we are reading from in the prognostic field
y_source_indexThe Y source index of the dimension we are reading from in the prognostic field
pid_locationLocation of the neighbouring process in the local stored data structures
current_pageThe current (next) buffer page to copy into
source_dataOptional source data which is read from

Definition at line 738 of file smagorinsky.F90.

740  type(model_state_type), intent(inout) :: current_state
741  integer, intent(in) :: corner_loc, pid_location, x_source_index, y_source_index
742  integer, intent(inout) :: current_page(:)
743  type(neighbour_description_type), intent(inout) :: neighbour_description
744  type(field_data_wrapper_type), dimension(:), intent(in), optional :: source_data
745 
746  call copy_corner_to_buffer(current_state%local_grid, neighbour_description%send_corner_buffer, &
747  current_state%diff_coefficient%data, corner_loc, x_source_index, y_source_index, current_page(pid_location))
748 
749  current_page(pid_location)=current_page(pid_location)+1
Here is the caller graph for this function:

◆ copy_diff_to_halo_buffer()

subroutine smagorinsky_mod::copy_diff_to_halo_buffer ( type(model_state_type), intent(inout)  current_state,
type(neighbour_description_type), intent(inout)  neighbour_description,
integer, intent(in)  dim,
integer, intent(in)  source_index,
integer, intent(in)  pid_location,
integer, dimension(:), intent(inout)  current_page,
type(field_data_wrapper_type), dimension(:), intent(in), optional  source_data 
)
private

Copies the diffusion field data to halo buffers for a specific process in a dimension and halo cell.

Parameters
current_stateThe current model state
neighbour_descriptionsDescription of the neighbour halo swapping status
dimDimension to copy from
source_indexThe source index of the dimension we are reading from in the prognostic field
pid_locationLocation of the neighbouring process in the local stored data structures
current_pageThe current (next) buffer page to copy into
source_dataOptional source data which is read from

Definition at line 715 of file smagorinsky.F90.

717  type(model_state_type), intent(inout) :: current_state
718  integer, intent(in) :: dim, pid_location, source_index
719  integer, intent(inout) :: current_page(:)
720  type(neighbour_description_type), intent(inout) :: neighbour_description
721  type(field_data_wrapper_type), dimension(:), intent(in), optional :: source_data
722 
723  call copy_field_to_buffer(current_state%local_grid, neighbour_description%send_halo_buffer, &
724  current_state%diff_coefficient%data, dim, source_index, current_page(pid_location))
725 
726  current_page(pid_location)=current_page(pid_location)+1
Here is the caller graph for this function:

◆ copy_vis_corners_to_halo_buffer()

subroutine smagorinsky_mod::copy_vis_corners_to_halo_buffer ( type(model_state_type), intent(inout)  current_state,
type(neighbour_description_type), intent(inout)  neighbour_description,
integer, intent(in)  corner_loc,
integer, intent(in)  x_source_index,
integer, intent(in)  y_source_index,
integer, intent(in)  pid_location,
integer, dimension(:), intent(inout)  current_page,
type(field_data_wrapper_type), dimension(:), intent(in), optional  source_data 
)
private

Copies the viscosity field corner data to halo buffers for a specific process.

Parameters
current_stateThe current model state
neighbour_descriptionsDescription of the neighbour halo swapping status
corner_locThe corner location
x_source_indexThe X source index of the dimension we are reading from in the prognostic field
y_source_indexThe Y source index of the dimension we are reading from in the prognostic field
pid_locationLocation of the neighbouring process in the local stored data structures
current_pageThe current (next) buffer page to copy into
source_dataOptional source data which is read from

Definition at line 693 of file smagorinsky.F90.

695  type(model_state_type), intent(inout) :: current_state
696  integer, intent(in) :: corner_loc, pid_location, x_source_index, y_source_index
697  integer, intent(inout) :: current_page(:)
698  type(neighbour_description_type), intent(inout) :: neighbour_description
699  type(field_data_wrapper_type), dimension(:), intent(in), optional :: source_data
700 
701  call copy_corner_to_buffer(current_state%local_grid, neighbour_description%send_corner_buffer, &
702  current_state%vis_coefficient%data, corner_loc, x_source_index, y_source_index, current_page(pid_location))
703 
704  current_page(pid_location)=current_page(pid_location)+1
Here is the caller graph for this function:

◆ copy_vis_to_halo_buffer()

subroutine smagorinsky_mod::copy_vis_to_halo_buffer ( type(model_state_type), intent(inout)  current_state,
type(neighbour_description_type), intent(inout)  neighbour_description,
integer, intent(in)  dim,
integer, intent(in)  source_index,
integer, intent(in)  pid_location,
integer, dimension(:), intent(inout)  current_page,
type(field_data_wrapper_type), dimension(:), intent(in), optional  source_data 
)
private

Copies the viscosity field data to halo buffers for a specific process in a dimension and halo cell.

Parameters
current_stateThe current model state
neighbour_descriptionsDescription of the neighbour halo swapping status
dimDimension to copy from
source_indexThe source index of the dimension we are reading from in the prognostic field
pid_locationLocation of the neighbouring process in the local stored data structures
current_pageThe current (next) buffer page to copy into
source_dataOptional source data which is read from

Definition at line 670 of file smagorinsky.F90.

672  type(model_state_type), intent(inout) :: current_state
673  integer, intent(in) :: dim, pid_location, source_index
674  integer, intent(inout) :: current_page(:)
675  type(neighbour_description_type), intent(inout) :: neighbour_description
676  type(field_data_wrapper_type), dimension(:), intent(in), optional :: source_data
677 
678  call copy_field_to_buffer(current_state%local_grid, neighbour_description%send_halo_buffer, &
679  current_state%vis_coefficient%data, dim, source_index, current_page(pid_location))
680 
681  current_page(pid_location)=current_page(pid_location)+1
Here is the caller graph for this function:

◆ finalisation_callback()

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

Called when the model is finishing up, will finalise the halo communications represented by the state.

Parameters
current_stateThe current model state

Definition at line 134 of file smagorinsky.F90.

135  type(model_state_type), target, intent(inout) :: current_state
136 
137  call finalise_halo_communication(current_state%viscosity_halo_swap_state)
138  call finalise_halo_communication(current_state%diffusion_halo_swap_state)
Here is the caller graph for this function:

◆ initialisation_callback()

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

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

Parameters
current_stateThe current model state

Definition at line 44 of file smagorinsky.F90.

45  type(model_state_type), target, intent(inout) :: current_state
46 
47 
48  if (.not. is_component_enabled(current_state%options_database, "diffusion")) then
49  call log_master_log(log_error, "Smagorinsky requires the diffusion component to be enabled")
50  end if
51  if (.not. is_component_enabled(current_state%options_database, "viscosity")) then
52  call log_master_log(log_error, "Smagorinsky requires the viscosity component to be enabled")
53  end if
54  if (.not. current_state%use_viscosity_and_diffusion) then
55  call log_master_log(log_error, "Smagorinsky requires use_viscosity_and_diffusion=.true. or monc will fail")
56  endif
57 
58  if (.not. is_component_enabled(current_state%options_database, "lower_bc")) then
59  call log_master_log(log_info, "LOWERBC is disabled, zero diff and vis_coeff on level 1")
60  current_state%vis_coefficient%data(1,:,:)=0.0_default_precision
61  current_state%diff_coefficient%data(1,:,:)=0.0_default_precision
62  endif
63 
64  eps=0.01_default_precision
65  repsh=0.5_default_precision/eps
66  thcona=ratio_mol_wts*current_state%thref0
67  thconb=thcona-current_state%thref0
68  thconap1=thcona
69 
70  subb=options_get_real(current_state%options_database, "smag-subb")
71  subc=options_get_real(current_state%options_database, "smag-subc")
72 
73  subg=1.2_default_precision
74  subh=0.0_default_precision
75  subr=4.0_default_precision
76  pr_n=0.7_default_precision
77  suba=1.0_default_precision/pr_n
78  ric=0.25_default_precision
79  ricinv=1.0_default_precision/ric
80  if (current_state%use_viscosity_and_diffusion) then
81  call init_halo_communication(current_state, get_single_field_per_halo_cell, &
82  current_state%viscosity_halo_swap_state, 2, .true.)
83  call init_halo_communication(current_state, get_single_field_per_halo_cell, &
84  current_state%diffusion_halo_swap_state, 2, .true.)
85  end if
86 
87  if (.not. current_state%passive_q) then
88  iqv = get_q_index(standard_q_names%VAPOUR, 'smagorinsky')
89  current_state%water_vapour_mixing_ratio_index=iqv
90 
91  iql = get_q_index(standard_q_names%CLOUD_LIQUID_MASS, 'smagorinsky')
92  current_state%liquid_water_mixing_ratio_index=iql
93  endif
94 
Here is the caller graph for this function:

◆ moist_ri_1()

real(kind=default_precision) function, dimension(current_state%local_grid%size(z_index)) smagorinsky_mod::moist_ri_1 ( type(model_state_type), intent(inout), target  current_state,
real(kind=default_precision), dimension(:), intent(in)  ssq,
type(prognostic_field_type), intent(inout)  th,
type(prognostic_field_type), dimension(:), intent(inout)  q 
)
private

Calculates numerator of Richarson number as the difference in buoyancy between a parcel at K+1 and a parcel lifted from K to K+1. Condensation is calculated by assuming that TL=T-(Lv/CP)*QL+gz/CP and QT = QV+QL are conserved on lifting. Microphysics conversions are assumed to happen slowly compared to turbulent processes so are, therefore, neglected in the calculation.

Parameters
current_stateThe current model state
ssqThe half squared strain rate
Returns
The Richardson number calculated here for a specific column

Definition at line 334 of file smagorinsky.F90.

335  type(model_state_type), target, intent(inout) :: current_state
336  real(kind=default_precision), dimension(:), intent(in) :: ssq
337  type(prognostic_field_type), intent(inout) :: th
338  type(prognostic_field_type), dimension(:), intent(inout) :: q
339  real(kind=default_precision), dimension(current_state%local_grid%size(Z_INDEX)) :: moist_ri_1
340 
341  integer :: k
342  real(kind=default_precision) :: liquid_water_at_kp1, temperature_at_kp1, dtlref, tmpinv, &
343  liquid_water_static_energy_temperature, total_water_substance, total_water_substance_p1, scltmp, qloading
344 
345  ! Set up the liquid water static energy & total water
346  do k=2, current_state%local_grid%size(z_index)
347  liquid_water_static_energy_temperature = &
348  th%data(k, current_state%column_local_y, current_state%column_local_x)*&
349  current_state%global_grid%configuration%vertical%rprefrcp(k) - rlvap_over_cp*&
350  q(current_state%liquid_water_mixing_ratio_index)%data(k, current_state%column_local_y, current_state%column_local_x)
351  if (k==2) then
352  total_water_substance = q(current_state%water_vapour_mixing_ratio_index)%data(&
353  k, current_state%column_local_y, current_state%column_local_x) + &
354  q(current_state%liquid_water_mixing_ratio_index)%data(k, current_state%column_local_y, current_state%column_local_x)
355  else
356  total_water_substance=total_water_substance_p1
357  end if
358  total_water_substance_p1 = q(current_state%water_vapour_mixing_ratio_index)%data(&
359  k+1, current_state%column_local_y, current_state%column_local_x) + &
360  q(current_state%liquid_water_mixing_ratio_index)%data(k+1, current_state%column_local_y, current_state%column_local_x)
361  if (k .le. current_state%local_grid%size(z_index)-1) then
362  tmpinv = 1.0_default_precision/ssq(k)
363  dtlref = current_state%global_grid%configuration%vertical%tref(k+1)-&
364  current_state%global_grid%configuration%vertical%tref(k)+(g/cp)*&
365  current_state%global_grid%configuration%vertical%dzn(k+1)
366  temperature_at_kp1=current_state%global_grid%configuration%vertical%tstarpr(k+1)+&
367  current_state%global_grid%configuration%vertical%tref(k+1)+&
368  (liquid_water_static_energy_temperature+rlvap_over_cp*(total_water_substance-&
369  current_state%global_grid%configuration%vertical%qsat(k+1)) - &
370  dtlref - current_state%global_grid%configuration%vertical%tstarpr(k+1) )*&
371  current_state%global_grid%configuration%vertical%qsatfac(k+1)
372  liquid_water_at_kp1=total_water_substance -( current_state%global_grid%configuration%vertical%qsat(k+1)+&
373  current_state%global_grid%configuration%vertical%dqsatdt(k+1)*(temperature_at_kp1- &
374  current_state%global_grid%configuration%vertical%tstarpr(k+1)-&
375  current_state%global_grid%configuration%vertical%tref(k+1)) )
376  if (liquid_water_at_kp1 .le. 0.0_default_precision) then
377  liquid_water_at_kp1=0.0_default_precision
378  temperature_at_kp1=current_state%global_grid%configuration%vertical%tref(k+1)+&
379  (liquid_water_static_energy_temperature-dtlref)
380  end if
381  scltmp=current_state%global_grid%configuration%vertical%rdzn(k+1)*&
382  current_state%global_grid%configuration%vertical%buoy_co(k)
383  ! Calculate the Richardson number taking into account the water
384  ! loading terms from the other hydrometeors in the buoyancy
385  ! calculation
386  qloading = current_state%cq(current_state%water_vapour_mixing_ratio_index)*&
387  (total_water_substance_p1-total_water_substance)-ratio_mol_wts*&
388  (q(current_state%liquid_water_mixing_ratio_index)%data(&
389  k+1, current_state%column_local_y, current_state%column_local_x)-liquid_water_at_kp1)
390  !if(IRAINP.gt.0) qloading = qloading + current_state%cq(IQR)*(l_q(J,K+1,IQR)-l_q(J,K,IQR))
391  !if(ISNOWP.gt.0) qloading = qloading + current_state%cq(IQS)*(l_q(J,K+1,IQS)-l_q(J,K,IQS))
392  !if(ICECLP.gt.0) qloading = qloading + current_state%cq(IQI)*(l_q(J,K+1,IQI)-l_q(J,K,IQI))
393  !if(IGRAUP.gt.0) qloading = qloading + currnet_state%cq(IQG)*(l_q(J,K+1,IQG)-l_q(J,K,IQG))
394 
395  moist_ri_1(k) = ((th%data(k+1, current_state%column_local_y, current_state%column_local_x)+&
396  current_state%global_grid%configuration%vertical%thref(k+1) - temperature_at_kp1*&
397  current_state%global_grid%configuration%vertical%prefrcp(k+1))&
398  +current_state%global_grid%configuration%vertical%thref(k+1)*qloading)*scltmp*tmpinv
399  end if
400  end do
Here is the caller graph for this function:

◆ moist_ri_2()

real(kind=default_precision) function, dimension(current_state%local_grid%size(z_index)) smagorinsky_mod::moist_ri_2 ( type(model_state_type), intent(inout), target  current_state,
real(kind=default_precision), dimension(:), intent(in)  ssq,
type(prognostic_field_type), intent(inout)  th,
type(prognostic_field_type), dimension(:), intent(inout)  q 
)
private

Calculates another "moist version" of the Richardson number based on "change in subgrid buoyancy flux when a fraction EPS is exchanged".

Parameters
current_stateThe current model state
ssqThe half squared strain rate
Returns
The Richardson number calculated here for a specific column

Definition at line 239 of file smagorinsky.F90.

240  type(model_state_type), target, intent(inout) :: current_state
241  real(kind=default_precision), dimension(:), intent(in) :: ssq
242  type(prognostic_field_type), intent(inout) :: th
243  type(prognostic_field_type), dimension(:), intent(inout) :: q
244  real(kind=default_precision), dimension(current_state%local_grid%size(Z_INDEX)) :: moist_ri_2
245 
246  integer :: k
247  real(kind=default_precision) :: tmpinv, scltmp, qt_k, qt_kp1, thlpr_k, thlpr_kp1, qlli, qlui, &
248  thvli, thvui, thlpr_un, thlpr_ln, qtun, qtln, thvln, thvun, qlln, qlun
249 
250  do k=2,current_state%local_grid%size(z_index)-1
251  tmpinv = 1.0_default_precision/ssq(k)
252  scltmp=current_state%global_grid%configuration%vertical%rdzn(k+1)*&
253  current_state%global_grid%configuration%vertical%buoy_co(k)
254  if (current_state%use_anelastic_equations) then
255  thcona=ratio_mol_wts*current_state%global_grid%configuration%vertical%thref(k)
256  thconb=0.5_default_precision*(ratio_mol_wts-1.0_default_precision)*&
257  (current_state%global_grid%configuration%vertical%thref(k)+&
258  current_state%global_grid%configuration%vertical%thref(k+1))
259  thconap1=ratio_mol_wts*current_state%global_grid%configuration%vertical%thref(k+1)
260  end if
261 
262  qt_k = q(current_state%water_vapour_mixing_ratio_index)%data(&
263  k, current_state%column_local_y, current_state%column_local_x)+&
264  q(current_state%liquid_water_mixing_ratio_index)%data(k, current_state%column_local_y, current_state%column_local_x)
265  qt_kp1 = q(current_state%water_vapour_mixing_ratio_index)%data(&
266  k+1, current_state%column_local_y, current_state%column_local_x)+&
267  q(current_state%liquid_water_mixing_ratio_index)%data(k+1, current_state%column_local_y, current_state%column_local_x)
268  thlpr_k = th%data(k, current_state%column_local_y, current_state%column_local_x) -&
269  rlvap_over_cp*q(current_state%liquid_water_mixing_ratio_index)%data(&
270  k, current_state%column_local_y, current_state%column_local_x)*&
271  current_state%global_grid%configuration%vertical%prefrcp(k)
272  ! calculate theta_l at constant (level K) pressure
273  thlpr_kp1 = th%data(k+1, current_state%column_local_y, current_state%column_local_x) - &
274  rlvap_over_cp*q(current_state%liquid_water_mixing_ratio_index)%data(&
275  k+1, current_state%column_local_y, current_state%column_local_x)*&
276  current_state%global_grid%configuration%vertical%prefrcp(k)
277 
278  !.......CALCULATE QL AND BUOYANCY OF LEVELS K AND K+1 (LOWER AND UPPER)
279  !.......QL CALC AT CONSTANT PRESSURE
280  qlli = max(0.0_default_precision, (qt_k-(current_state%global_grid%configuration%vertical%qsat(k) +&
281  current_state%global_grid%configuration%vertical%dqsatdt(k)*&
282  (current_state%global_grid%configuration%vertical%rprefrcp(k)*thlpr_k-&
283  current_state%global_grid%configuration%vertical%tstarpr(k))))*&
284  current_state%global_grid%configuration%vertical%qsatfac(k))
285  qlui = max(0.0_default_precision, (qt_kp1-(current_state%global_grid%configuration%vertical%qsat(k) +&
286  current_state%global_grid%configuration%vertical%dqsatdt(k)*&
287  (current_state%global_grid%configuration%vertical%rprefrcp(k)*(thlpr_kp1+&
288  current_state%global_grid%configuration%vertical%thref(k+1))-&
289  current_state%global_grid%configuration%vertical%tstarpr(k)-&
290  current_state%global_grid%configuration%vertical%tref(k))))*&
291  current_state%global_grid%configuration%vertical%qsatfac(k))
292  thvli = thlpr_k+current_state%global_grid%configuration%vertical%thref(k) +&
293  ( rlvap_over_cp*current_state%global_grid%configuration%vertical%prefrcp(k)-thcona ) * qlli
294  thvui = thlpr_kp1+current_state%global_grid%configuration%vertical%thref(k+1) +&
295  ( rlvap_over_cp*current_state%global_grid%configuration%vertical%prefrcp(k)-thconap1 ) * qlui
296 
297  !.......MIX CONSERVATIVE QUANTITIES AT CONSTANT PRESSURE
298  !.......NOTE THAT FOR COMPUTATIONAL CONVENIENCE THE RELEVANT LEVEL'S
299  !.......THREF HAS BEEN OMITTED FROM ALL THV'S, THLPR_LN AND THLPR_UN AS
300  !.......ONLY THE DIFFERENCE BETWEEN THV'S AT THE SAME LEVEL IS REQUIRED.
301  thlpr_un = thlpr_kp1 - eps*(thlpr_kp1 - thlpr_k + current_state%global_grid%configuration%vertical%dthref(k))
302  thlpr_ln = thlpr_k + eps*(thlpr_kp1 - thlpr_k + current_state%global_grid%configuration%vertical%dthref(k))
303  qtun = qt_kp1 - eps*(qt_kp1 - qt_k)
304  qtln = qt_k + eps*(qt_kp1 - qt_k)
305  !.......CALCULATE LIQUID WATER AND BUOYANCY OF MIXED AIR
306  qlln = max(0.0_default_precision,( qtln - ( current_state%global_grid%configuration%vertical%qsat(k) +&
307  current_state%global_grid%configuration%vertical%dqsatdt(k)*&
308  (current_state%global_grid%configuration%vertical%rprefrcp(k)*&
309  thlpr_ln-current_state%global_grid%configuration%vertical%tstarpr(k)))&
310  )*current_state%global_grid%configuration%vertical%qsatfac(k))
311  qlun = max(0.0_default_precision,( qtun - ( current_state%global_grid%configuration%vertical%qsat(k) + &
312  current_state%global_grid%configuration%vertical%dqsatdt(k)*&
313  (current_state%global_grid%configuration%vertical%rprefrcp(k)*&
314  (thlpr_un+current_state%global_grid%configuration%vertical%thref(k+1))-&
315  current_state%global_grid%configuration%vertical%tstarpr(k)-current_state%global_grid%configuration%vertical%tref(k)))&
316  )*current_state%global_grid%configuration%vertical%qsatfac(k))
317  thvln = thlpr_ln+current_state%global_grid%configuration%vertical%thref(k) +&
318  ( rlvap_over_cp*current_state%global_grid%configuration%vertical%prefrcp(k)-thcona) * qlln
319  thvun = thlpr_un+current_state%global_grid%configuration%vertical%thref(k+1) +&
320  ( rlvap_over_cp*current_state%global_grid%configuration%vertical%prefrcp(k)-thconap1) * qlun
321 
322  !.....CHANGE IN BUOYANCY DUE TO FRAC MIXING IS USED TO DETERMINE l_ri
323  moist_ri_2(k) = scltmp*tmpinv*(thconb*(qt_kp1-qt_k) +((thvln-thvli)-(thvun-thvui))*repsh)
324  end do
Here is the caller graph for this function:

◆ setfri()

subroutine smagorinsky_mod::setfri ( type(model_state_type), intent(inout), target  current_state,
real(kind=default_precision), dimension(:), intent(in)  richardson_number,
real(kind=default_precision), dimension(:), intent(in)  ssq 
)
private

Calculates the eddy viscosity (VIS) and diffusivity (DIFF) depending on the Richardson Number (RI) and half squared strain rate.

Parameters
current_stateThe current model state
richardson_numberRichardson number
ssqThe half squared strain rate

Definition at line 164 of file smagorinsky.F90.

165  type(model_state_type), target, intent(inout) :: current_state
166  real(kind=default_precision), dimension(:), intent(in) :: richardson_number, ssq
167 
168  integer :: k
169  real(kind=default_precision) :: rifac, sctmp
170 
171  do k=2,current_state%local_grid%size(z_index)-1
172  if (richardson_number(k) .le. 0.0_default_precision) then
173  current_state%vis_coefficient%data(k, current_state%column_local_y, current_state%column_local_x)=&
174  sqrt(1.0_default_precision-subc*richardson_number(k))
175  current_state%diff_coefficient%data(k, current_state%column_local_y, current_state%column_local_x)=&
176  suba*sqrt(1.-subb*richardson_number(k))
177  else if ((richardson_number(k) .gt. 0.0_default_precision) .and. (richardson_number(k) .lt. ric)) then
178  rifac=(1.-richardson_number(k)*ricinv)**4
179  current_state%vis_coefficient%data(k, current_state%column_local_y, current_state%column_local_x)=&
180  rifac*(1.-subh*richardson_number(k))
181  current_state%diff_coefficient%data(k, current_state%column_local_y, current_state%column_local_x)=&
182  rifac*suba*(1.-subg*richardson_number(k))
183  else
184  current_state%vis_coefficient%data(k, current_state%column_local_y, current_state%column_local_x)=0.0_default_precision
185  current_state%diff_coefficient%data(k, current_state%column_local_y, current_state%column_local_x)=0.0_default_precision
186  end if
187  sctmp=current_state%global_grid%configuration%vertical%rneutml_sq(k)*sqrt(ssq(k))
188  current_state%vis_coefficient%data(k, current_state%column_local_y, current_state%column_local_x)=&
189  current_state%vis_coefficient%data(k, current_state%column_local_y, current_state%column_local_x)*sctmp
190  current_state%diff_coefficient%data(k, current_state%column_local_y, current_state%column_local_x)=&
191  current_state%diff_coefficient%data(k, current_state%column_local_y, current_state%column_local_x)*sctmp
192  end do
193 
194  current_state%vis_coefficient%data(current_state%local_grid%size(z_index), &
195  current_state%column_local_y, current_state%column_local_x)= 0.0_default_precision
196  current_state%diff_coefficient%data(current_state%local_grid%size(z_index), &
197  current_state%column_local_y, current_state%column_local_x)= 0.0_default_precision
Here is the caller graph for this function:

◆ smagorinsky_get_descriptor()

type(component_descriptor_type) function, public smagorinsky_mod::smagorinsky_get_descriptor

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

Returns
The termination check component descriptor

Definition at line 34 of file smagorinsky.F90.

35  smagorinsky_get_descriptor%name="smagorinsky"
36  smagorinsky_get_descriptor%version=0.1
37  smagorinsky_get_descriptor%initialisation=>initialisation_callback
38  smagorinsky_get_descriptor%timestep=>timestep_callback
39  smagorinsky_get_descriptor%finalisation=>finalisation_callback
Here is the call graph for this function:

◆ timestep_callback()

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

For each none halo cell this will calculate the subgrid terms for viscosity and diffusion.

Parameters
current_stateThe current model state

Definition at line 99 of file smagorinsky.F90.

100  type(model_state_type), target, intent(inout) :: current_state
101 
102  real(kind=default_precision), dimension(current_state%local_grid%size(Z_INDEX)) :: richardson_number, ssq
103 
104  if (.not. current_state%halo_column) then
105  if (.not. current_state%use_viscosity_and_diffusion) then
106  current_state%vis_coefficient%data(2:,:,:)=0.0_default_precision
107  current_state%diff_coefficient%data(2:,:,:)=0.0_default_precision
108  current_state%cvis=0.0_default_precision
109  else
110  if (current_state%field_stepping == forward_stepping) then
111  ssq=calculate_half_squared_strain_rate(current_state, current_state%u, current_state%v, current_state%w)
112  richardson_number=calculate_richardson_number(current_state, ssq, current_state%th, current_state%q)
113  else
114  ssq=calculate_half_squared_strain_rate(current_state, current_state%zu, current_state%zv, current_state%zw)
115  richardson_number=calculate_richardson_number(current_state, ssq, current_state%zth, current_state%zq)
116  end if
117  call setfri(current_state, richardson_number, ssq)
118  if (is_component_enabled(current_state%options_database, "cfltest")) then
119  call update_viscous_number(current_state)
120  endif
121  end if
122  end if
123  if (current_state%last_timestep_column) then
124  ! need to ensure not already in progress
125  call initiate_nonblocking_halo_swap(current_state, current_state%diffusion_halo_swap_state, &
126  copy_diff_to_halo_buffer, copy_diff_corners_to_halo_buffer)
127  call initiate_nonblocking_halo_swap(current_state, current_state%viscosity_halo_swap_state, &
128  copy_vis_to_halo_buffer, copy_vis_corners_to_halo_buffer)
129  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ update_viscous_number()

subroutine smagorinsky_mod::update_viscous_number ( type(model_state_type), intent(inout), target  current_state)
private

Update viscous number based upon the viscosity and diffusivity coefficients.

Parameters
current_stateThe current model state

Definition at line 143 of file smagorinsky.F90.

144  type(model_state_type), target, intent(inout) :: current_state
145 
146  integer :: k
147 
148  if (mod(current_state%timestep, current_state%cfl_frequency) == 1 .or. &
149  current_state%timestep-current_state%start_timestep .le. current_state%cfl_frequency) then
150  do k=2, current_state%local_grid%size(z_index)-1
151  current_state%cvis=max(current_state%cvis, max(current_state%vis_coefficient%data(k, current_state%column_local_y, &
152  current_state%column_local_x),current_state%diff_coefficient%data(k, current_state%column_local_y, &
153  current_state%column_local_x))*(current_state%global_grid%configuration%vertical%rdzn(k+1)**2+&
154  current_state%global_grid%configuration%horizontal%cx2+current_state%global_grid%configuration%horizontal%cy2))
155  end do
156  end if
Here is the caller graph for this function:

Variable Documentation

◆ eps

real(kind=default_precision) smagorinsky_mod::eps
private

Definition at line 23 of file smagorinsky.F90.

23  real(kind=default_precision) :: eps, repsh, thcona, thconb, thconap1, suba, subb, subc, subg, subh, subr, pr_n, ric, ricinv

◆ iql

integer smagorinsky_mod::iql
private

Definition at line 28 of file smagorinsky.F90.

◆ iqv

integer smagorinsky_mod::iqv
private

Definition at line 28 of file smagorinsky.F90.

28  integer :: iqv, iql ! index for water vapour and liquid

◆ pr_n

real(kind=default_precision) smagorinsky_mod::pr_n
private

Definition at line 23 of file smagorinsky.F90.

◆ repsh

real(kind=default_precision) smagorinsky_mod::repsh
private

Definition at line 23 of file smagorinsky.F90.

◆ ric

real(kind=default_precision) smagorinsky_mod::ric
private

Definition at line 23 of file smagorinsky.F90.

◆ richardson_number_calculation

integer, parameter smagorinsky_mod::richardson_number_calculation =2
private

Definition at line 22 of file smagorinsky.F90.

22  integer, parameter :: RICHARDSON_NUMBER_CALCULATION=2

◆ ricinv

real(kind=default_precision) smagorinsky_mod::ricinv
private

Definition at line 23 of file smagorinsky.F90.

◆ suba

real(kind=default_precision) smagorinsky_mod::suba
private

Definition at line 23 of file smagorinsky.F90.

◆ subb

real(kind=default_precision) smagorinsky_mod::subb
private

Definition at line 23 of file smagorinsky.F90.

◆ subc

real(kind=default_precision) smagorinsky_mod::subc
private

Definition at line 23 of file smagorinsky.F90.

◆ subg

real(kind=default_precision) smagorinsky_mod::subg
private

Definition at line 23 of file smagorinsky.F90.

◆ subh

real(kind=default_precision) smagorinsky_mod::subh
private

Definition at line 23 of file smagorinsky.F90.

◆ subr

real(kind=default_precision) smagorinsky_mod::subr
private

Definition at line 23 of file smagorinsky.F90.

◆ thcona

real(kind=default_precision) smagorinsky_mod::thcona
private

Definition at line 23 of file smagorinsky.F90.

◆ thconap1

real(kind=default_precision) smagorinsky_mod::thconap1
private

Definition at line 23 of file smagorinsky.F90.

◆ thconb

real(kind=default_precision) smagorinsky_mod::thconb
private

Definition at line 23 of file smagorinsky.F90.

logging_mod::log_error
integer, parameter, public log_error
Only log ERROR messages.
Definition: logging.F90:11
logging_mod::log_info
integer, parameter, public log_info
Log INFO, WARNING and ERROR messages.
Definition: logging.F90:13
science_constants_mod::cp
real(kind=default_precision), public cp
Definition: scienceconstants.F90:13
grids_mod::z_index
integer, parameter, public z_index
Grid index parameters.
Definition: grids.F90:14
logging_mod::log_master_log
subroutine, public log_master_log(level, message)
Will log just from the master process.
Definition: logging.F90:47
optionsdatabase_mod::options_get_real
real(kind=default_precision) function, public options_get_real(options_database, key, index)
Retrieves a real value from the database that matches the provided key.
Definition: optionsdatabase.F90:91
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