MONC
Functions/Subroutines | Variables
viscosity_mod Module Reference

Computes the viscosity dynamics for the U,V,W source terms. More...

Functions/Subroutines

type(component_descriptor_type) function, public viscosity_get_descriptor ()
 Provides the descriptor back to the caller and is used in component registration. More...
 
subroutine field_information_retrieval_callback (current_state, name, field_information)
 Field information retrieval callback, this returns information for a specific components published field. More...
 
subroutine field_value_retrieval_callback (current_state, name, field_value)
 Field value retrieval callback, this returns the value of a specific published field. More...
 
subroutine initialisation_callback (current_state)
 Sets up the stencil_mod (used in interpolation) and allocates data for the flux fields. More...
 
subroutine finalisation_callback (current_state)
 
subroutine timestep_callback (current_state)
 At each timestep will compute the viscosity U,V,W source terms. More...
 
subroutine calculate_viscous_sources (current_state, local_y, local_x, tau12, tau12_ym1, tau12m1, tau11, tau22, tau22_yp1, tau33, tau11p1, tau13, tau13m1, tau23, tau23_ym1)
 Calculates the viscous sources based upon TAU for the U,V and W fields. More...
 
subroutine calculate_tau (current_state, local_y, local_x, zu, zv, zw, tau12, tau12_ym1, tau12m1, tau11, tau22, tau22_yp1, tau33, tau11p1, tau13, tau13m1, tau23, tau23_ym1)
 Calculated TAU which is used in computing the viscous source terms. More...
 
subroutine perform_local_data_copy_for_vis (current_state, halo_depth, involve_corners, source_data)
 Does local data copying for viscosity coefficient variable halo swap. More...
 
subroutine copy_halo_buffer_to_vis (current_state, neighbour_description, dim, target_index, neighbour_location, current_page, source_data)
 Copies the halo buffer to halo location for the viscosity coefficient field. More...
 
subroutine copy_halo_buffer_to_vis_corners (current_state, neighbour_description, corner_loc, x_target_index, y_target_index, neighbour_location, current_page, source_data)
 Copies the corner halo buffer to the viscosity coefficient field corners. More...
 
subroutine save_precomponent_tendencies (current_state, cxn, cyn, txn, tyn)
 Save the 3d tendencies coming into the component. More...
 
subroutine compute_component_tendencies (current_state, cxn, cyn, txn, tyn)
 Computation of component tendencies. More...
 
subroutine set_published_field_value (field_value, real_1d_field, real_2d_field, real_3d_field)
 Sets the published field value from the temporary diagnostic values held by this component. More...
 

Variables

real(kind=default_precision), dimension(:), allocatable u_viscosity
 
real(kind=default_precision), dimension(:), allocatable v_viscosity
 
real(kind=default_precision), dimension(:), allocatable w_viscosity
 
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_u
 
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_v
 
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_w
 
logical l_tend_3d_u
 
logical l_tend_3d_v
 
logical l_tend_3d_w
 
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_u
 
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_v
 
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_w
 
logical l_tend_pr_tot_u
 
logical l_tend_pr_tot_v
 
logical l_tend_pr_tot_w
 
integer diagnostic_generation_frequency
 

Detailed Description

Computes the viscosity dynamics for the U,V,W source terms.

Function/Subroutine Documentation

◆ calculate_tau()

subroutine viscosity_mod::calculate_tau ( type(model_state_type), intent(inout), target  current_state,
integer, intent(in)  local_y,
integer, intent(in)  local_x,
type(prognostic_field_type), intent(inout)  zu,
type(prognostic_field_type), intent(inout)  zv,
type(prognostic_field_type), intent(inout)  zw,
real(kind=default_precision), dimension(:), intent(out)  tau12,
real(kind=default_precision), dimension(:), intent(out)  tau12_ym1,
real(kind=default_precision), dimension(:), intent(out)  tau12m1,
real(kind=default_precision), dimension(:), intent(out)  tau11,
real(kind=default_precision), dimension(:), intent(out)  tau22,
real(kind=default_precision), dimension(:), intent(out)  tau22_yp1,
real(kind=default_precision), dimension(:), intent(out)  tau33,
real(kind=default_precision), dimension(:), intent(out)  tau11p1,
real(kind=default_precision), dimension(:), intent(out)  tau13,
real(kind=default_precision), dimension(:), intent(out)  tau13m1,
real(kind=default_precision), dimension(:), intent(out)  tau23,
real(kind=default_precision), dimension(:), intent(out)  tau23_ym1 
)
private

Calculated TAU which is used in computing the viscous source terms.

Parameters
current_stateThe current model state
local_yLocal Y index
local_xLocal X index

Definition at line 336 of file viscosity.F90.

338  type(model_state_type), target, intent(inout) :: current_state
339  type(prognostic_field_type), intent(inout) :: zu, zv, zw
340  integer, intent(in) :: local_y, local_x
341  real(kind=default_precision), dimension(:), intent(out) :: tau12, tau12m1, tau11, tau22, tau22_yp1, tau33, &
342  tau11p1, tau13, tau13m1, tau23, tau23_ym1, tau12_ym1
343 
344  integer :: k
345  real(kind=default_precision) :: vistmp, vis12, vis12a, visonp2, visonp2a, vis13, vis23
346 
347  ! Do p levels and w-levels
348  do k=1, current_state%local_grid%size(z_index)
349  if (k .gt. 1) then
350  vistmp=current_state%vis_coefficient%data(k, local_y, local_x)+current_state%vis_coefficient%data(k, local_y+1, local_x)+&
351  current_state%vis_coefficient%data(k-1, local_y, local_x)+current_state%vis_coefficient%data(k-1, local_y+1, local_x)
352  vis12=0.125_default_precision*(vistmp + current_state%vis_coefficient%data(k, local_y, local_x+1)+&
353  current_state%vis_coefficient%data(k, local_y+1, local_x+1)+&
354  current_state%vis_coefficient%data(k-1, local_y, local_x+1)+&
355  current_state%vis_coefficient%data(k-1, local_y+1, local_x+1))
356  tau12(k)=0.0_default_precision
357 #ifdef U_ACTIVE
358  tau12(k)=(zu%data(k, local_y+1, local_x)-zu%data(k, local_y, local_x))*&
359  current_state%global_grid%configuration%horizontal%cy
360 #endif
361 #ifdef V_ACTIVE
362  tau12(k)=tau12(k)+(zv%data(k, local_y, local_x+1)-zv%data(k, local_y, local_x))*&
363  current_state%global_grid%configuration%horizontal%cx
364 #endif
365  tau12(k)=tau12(k)*current_state%global_grid%configuration%vertical%rhon(k)*vis12
366 
367  vis12a=0.125_default_precision*(vistmp +current_state%vis_coefficient%data(k, local_y, local_x-1)+&
368  current_state%vis_coefficient%data(k, local_y+1, local_x-1)+&
369  current_state%vis_coefficient%data(k-1, local_y, local_x-1)+&
370  current_state%vis_coefficient%data(k-1, local_y+1, local_x-1))
371  tau12m1(k)=0.0_default_precision
372 #ifdef U_ACTIVE
373  tau12m1(k)=(zu%data(k, local_y+1, local_x-1)-zu%data(k, local_y, local_x-1))*&
374  current_state%global_grid%configuration%horizontal%cy
375 #endif
376 #ifdef V_ACTIVE
377  tau12m1(k)=tau12m1(k)+(zv%data(k, local_y, local_x)-zv%data(k, local_y, local_x-1))*&
378  current_state%global_grid%configuration%horizontal%cx
379 #endif
380  tau12m1(k)=tau12m1(k)*current_state%global_grid%configuration%vertical%rhon(k)*vis12a
381 
382  vistmp=current_state%vis_coefficient%data(k, local_y-1, local_x)+current_state%vis_coefficient%data(k, local_y, local_x)+&
383  current_state%vis_coefficient%data(k-1, local_y-1, local_x)+current_state%vis_coefficient%data(k-1, local_y, local_x)
384  vis12=0.125_default_precision*(vistmp + current_state%vis_coefficient%data(k, local_y-1, local_x+1)+&
385  current_state%vis_coefficient%data(k, local_y, local_x+1)+&
386  current_state%vis_coefficient%data(k-1, local_y-1, local_x+1)+&
387  current_state%vis_coefficient%data(k-1, local_y, local_x+1))
388  tau12_ym1(k)=0.0_default_precision
389 #ifdef U_ACTIVE
390  tau12_ym1(k)=(zu%data(k, local_y, local_x)-zu%data(k, local_y-1, local_x))*&
391  current_state%global_grid%configuration%horizontal%cy
392 #endif
393 #ifdef V_ACTIVE
394  tau12_ym1(k)=tau12_ym1(k)+(zv%data(k, local_y-1, local_x+1)-zv%data(k, local_y-1, local_x))*&
395  current_state%global_grid%configuration%horizontal%cx
396 #endif
397  tau12_ym1(k)=tau12_ym1(k)*current_state%global_grid%configuration%vertical%rhon(k)*vis12
398 
399  visonp2=current_state%global_grid%configuration%vertical%rhon(k)*&
400  (current_state%vis_coefficient%data(k, local_y, local_x)+current_state%vis_coefficient%data(k-1, local_y, local_x))
401 #ifdef U_ACTIVE
402  tau11(k)=visonp2*(zu%data(k, local_y, local_x)-zu%data(k, local_y, local_x-1))*&
403  current_state%global_grid%configuration%horizontal%cx
404 #else
405  tau11(k)=0.0_default_precision
406 #endif
407 #ifdef V_ACTIVE
408  tau22(k)=visonp2*(zv%data(k, local_y, local_x)-zv%data(k, local_y-1, local_x))*&
409  current_state%global_grid%configuration%horizontal%cy
410 #else
411  tau22(k)=0.0_default_precision
412 #endif
413 #ifdef W_ACTIVE
414  tau33(k)=visonp2*(zw%data(k, local_y, local_x)-zw%data(k-1, local_y, local_x))*&
415  current_state%global_grid%configuration%vertical%rdz(k)
416 #else
417  tau33(k)=0.0_default_precision
418 #endif
419 
420 #ifdef V_ACTIVE
421  visonp2=current_state%global_grid%configuration%vertical%rhon(k)*&
422  (current_state%vis_coefficient%data(k, local_y+1, local_x)+&
423  current_state%vis_coefficient%data(k-1, local_y+1, local_x))
424  tau22_yp1(k)=visonp2*(zv%data(k, local_y+1, local_x)-zv%data(k, local_y, local_x))*&
425  current_state%global_grid%configuration%horizontal%cy
426 #endif
427 
428 #ifdef U_ACTIVE
429  visonp2a=current_state%global_grid%configuration%vertical%rhon(k)*&
430  (current_state%vis_coefficient%data(k, local_y, local_x+1)+&
431  current_state%vis_coefficient%data(k-1, local_y, local_x+1))
432  tau11p1(k)=visonp2a*(current_state%zu%data(k, local_y, local_x+1)-&
433  zu%data(k, local_y, local_x))*current_state%global_grid%configuration%horizontal%cx
434 #endif
435  else
436  tau12(k)=0.0_default_precision
437  tau12_ym1(k)=0.0_default_precision
438  tau12m1(k)=0.0_default_precision
439  tau11(k)=0.0_default_precision
440  tau22(k)=0.0_default_precision
441  tau22_yp1(k)=0.0_default_precision
442  tau33(k)=0.0_default_precision
443  tau11p1(k)=0.0_default_precision
444  end if
445  if (k .lt. current_state%local_grid%size(z_index)) then
446  vis13=0.5_default_precision*(current_state%vis_coefficient%data(k, local_y, local_x)+&
447  current_state%vis_coefficient%data(k, local_y, local_x+1))
448  tau13(k)=0.0_default_precision
449 #ifdef U_ACTIVE
450  tau13(k)=(zu%data(k+1, local_y, local_x)-zu%data(k, local_y, local_x))*&
451  current_state%global_grid%configuration%vertical%rdzn(k+1)
452 #endif
453 #ifdef W_ACTIVE
454  tau13(k)=tau13(k)+(zw%data(k, local_y, local_x+1)-zw%data(k, local_y, local_x))*&
455  current_state%global_grid%configuration%horizontal%cx
456 #endif
457  tau13(k)=tau13(k)*current_state%global_grid%configuration%vertical%rho(k)*vis13
458 
459  vis13=0.5_default_precision*(current_state%vis_coefficient%data(k, local_y, local_x-1)+&
460  current_state%vis_coefficient%data(k, local_y, local_x))
461  tau13m1(k)=0.0_default_precision
462 #ifdef U_ACTIVE
463  tau13m1(k)=(zu%data(k+1, local_y, local_x-1)-zu%data(k, local_y, local_x-1))*&
464  current_state%global_grid%configuration%vertical%rdzn(k+1)
465 #endif
466 #ifdef W_ACTIVE
467  tau13m1(k)=tau13m1(k)+(zw%data(k, local_y, local_x)-zw%data(k, local_y, local_x-1))*&
468  current_state%global_grid%configuration%horizontal%cx
469 #endif
470  tau13m1(k)=tau13m1(k)*current_state%global_grid%configuration%vertical%rho(k)*vis13
471 
472  vis23=0.5_default_precision*(current_state%vis_coefficient%data(k, local_y, local_x)+&
473  current_state%vis_coefficient%data(k, local_y+1, local_x))
474  tau23(k)=0.0_default_precision
475 #ifdef W_ACTIVE
476  tau23(k)=(zw%data(k, local_y+1, local_x)-zw%data(k, local_y, local_x))*&
477  current_state%global_grid%configuration%horizontal%cy
478 #endif
479 #ifdef V_ACTIVE
480  tau23(k)=tau23(k)+(zv%data(k+1, local_y, local_x)-zv%data(k, local_y, local_x))*&
481  current_state%global_grid%configuration%vertical%rdzn(k+1)
482 #endif
483  tau23(k)=tau23(k)*current_state%global_grid%configuration%vertical%rho(k)*vis23
484 
485  vis23=0.5_default_precision*(current_state%vis_coefficient%data(k, local_y-1, local_x)+&
486  current_state%vis_coefficient%data(k, local_y, local_x))
487  tau23_ym1(k)=0.0_default_precision
488 #ifdef W_ACTIVE
489  tau23_ym1(k)=(zw%data(k, local_y, local_x)-zw%data(k, local_y-1, local_x))*&
490  current_state%global_grid%configuration%horizontal%cy
491 #endif
492 #ifdef V_ACTIVE
493  tau23_ym1(k)=tau23_ym1(k)+(zv%data(k+1, local_y-1, local_x)-zv%data(k, local_y-1, local_x))*&
494  current_state%global_grid%configuration%vertical%rdzn(k+1)
495 #endif
496  tau23_ym1(k)=tau23_ym1(k)*current_state%global_grid%configuration%vertical%rho(k)*vis23
497  else
498  tau13(k)=0.0_default_precision
499  tau13m1(k)=0.0_default_precision
500  tau23(k)=0.0_default_precision
501  tau23_ym1(k)=0.0_default_precision
502  end if
503  end do
Here is the caller graph for this function:

◆ calculate_viscous_sources()

subroutine viscosity_mod::calculate_viscous_sources ( type(model_state_type), intent(inout), target  current_state,
integer, intent(in)  local_y,
integer, intent(in)  local_x,
real(kind=default_precision), dimension(:), intent(in)  tau12,
real(kind=default_precision), dimension(:), intent(in)  tau12_ym1,
real(kind=default_precision), dimension(:), intent(in)  tau12m1,
real(kind=default_precision), dimension(:), intent(in)  tau11,
real(kind=default_precision), dimension(:), intent(in)  tau22,
real(kind=default_precision), dimension(:), intent(in)  tau22_yp1,
real(kind=default_precision), dimension(:), intent(in)  tau33,
real(kind=default_precision), dimension(:), intent(in)  tau11p1,
real(kind=default_precision), dimension(:), intent(in)  tau13,
real(kind=default_precision), dimension(:), intent(in)  tau13m1,
real(kind=default_precision), dimension(:), intent(in)  tau23,
real(kind=default_precision), dimension(:), intent(in)  tau23_ym1 
)
private

Calculates the viscous sources based upon TAU for the U,V and W fields.

Parameters
current_stateThe current model state
local_yLocal Y index
local_xLocal X index

Definition at line 299 of file viscosity.F90.

301  type(model_state_type), target, intent(inout) :: current_state
302  integer, intent(in) :: local_y, local_x
303  real(kind=default_precision), dimension(:), intent(in) :: tau12, tau12_ym1, tau12m1, tau11, tau22, tau22_yp1, tau33, &
304  tau11p1, tau13, tau13m1, tau23, tau23_ym1
305 
306  integer :: k
307 
308  do k=2, current_state%local_grid%size(z_index)
309 #ifdef U_ACTIVE
310  u_viscosity(k)=((tau11p1(k)-tau11(k))*current_state%global_grid%configuration%horizontal%cx+(tau12(k)-tau12_ym1(k))*&
311  current_state%global_grid%configuration%horizontal%cy+(tau13(k)-tau13(k-1))*&
312  current_state%global_grid%configuration%vertical%rdz(k))/current_state%global_grid%configuration%vertical%rhon(k)
313  current_state%su%data(k, local_y, local_x)=current_state%su%data(k, local_y, local_x)+u_viscosity(k)
314 #endif
315 #ifdef V_ACTIVE
316  v_viscosity(k)=((tau12(k)-tau12m1(k))*current_state%global_grid%configuration%horizontal%cx+(tau22_yp1(k)-tau22(k))*&
317  current_state%global_grid%configuration%horizontal%cy+(tau23(k)-tau23(k-1))*&
318  current_state%global_grid%configuration%vertical%rdz(k))/current_state%global_grid%configuration%vertical%rhon(k)
319  current_state%sv%data(k, local_y, local_x)=current_state%sv%data(k, local_y, local_x)+v_viscosity(k)
320 #endif
321  end do
322 #ifdef W_ACTIVE
323  do k=2, current_state%local_grid%size(z_index)-1
324  w_viscosity(k)=((tau13(k)-tau13m1(k))*current_state%global_grid%configuration%horizontal%cx+(tau23(k)-tau23_ym1(k))*&
325  current_state%global_grid%configuration%horizontal%cy+(tau33(k+1)-tau33(k))*&
326  current_state%global_grid%configuration%vertical%rdzn(k+1))/current_state%global_grid%configuration%vertical%rho(k)
327  current_state%sw%data(k, local_y, local_x)=current_state%sw%data(k, local_y, local_x)+w_viscosity(k)
328  end do
329 #endif
Here is the caller graph for this function:

◆ compute_component_tendencies()

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

Computation of component tendencies.

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

Definition at line 594 of file viscosity.F90.

595  type(model_state_type), target, intent(inout) :: current_state
596  integer, intent(in) :: cxn, cyn, txn, tyn
597 
598  ! Calculate change in tendency due to component
599  if (l_tend_3d_u) then
600  tend_3d_u(:,tyn,txn)=current_state%su%data(:,cyn,cxn) - tend_3d_u(:,tyn,txn)
601  endif
602  if (l_tend_3d_v) then
603  tend_3d_v(:,tyn,txn)=current_state%sv%data(:,cyn,cxn) - tend_3d_v(:,tyn,txn)
604  endif
605  if (l_tend_3d_w) then
606  tend_3d_w(:,tyn,txn)=current_state%sw%data(:,cyn,cxn) - tend_3d_w(:,tyn,txn)
607  endif
608 
609  ! Add local tendency fields to the profile total
610  if (l_tend_pr_tot_u) then
611  tend_pr_tot_u(:)=tend_pr_tot_u(:) + tend_3d_u(:,tyn,txn)
612  endif
613  if (l_tend_pr_tot_v) then
614  tend_pr_tot_v(:)=tend_pr_tot_v(:) + tend_3d_v(:,tyn,txn)
615  endif
616  if (l_tend_pr_tot_w) then
617  tend_pr_tot_w(:)=tend_pr_tot_w(:) + tend_3d_w(:,tyn,txn)
618  endif
619 
Here is the caller graph for this function:

◆ copy_halo_buffer_to_vis()

subroutine viscosity_mod::copy_halo_buffer_to_vis ( type(model_state_type), intent(inout)  current_state,
type(neighbour_description_type), intent(inout)  neighbour_description,
integer, intent(in)  dim,
integer, intent(in)  target_index,
integer, intent(in)  neighbour_location,
integer, dimension(:), intent(inout)  current_page,
type(field_data_wrapper_type), dimension(:), intent(in), optional  source_data 
)
private

Copies the halo buffer to halo location for the viscosity coefficient field.

Parameters
current_stateThe current model state
neighbour_descriptionThe halo swapping description of the neighbour we are accessing the buffer of
dimThe dimension we receive for
target_indexThe target index for the dimension we are receiving for
neighbour_locationThe location in the local neighbour data stores of this neighbour
current_pageThe current, next, halo swap page to read from (all previous have been read and copied already)
source_dataOptional source data which is written into

Definition at line 527 of file viscosity.F90.

529  type(model_state_type), intent(inout) :: current_state
530  integer, intent(in) :: dim, target_index, neighbour_location
531  integer, intent(inout) :: current_page(:)
532  type(neighbour_description_type), intent(inout) :: neighbour_description
533  type(field_data_wrapper_type), dimension(:), intent(in), optional :: source_data
534 
535  call copy_buffer_to_field(current_state%local_grid, neighbour_description%recv_halo_buffer, &
536  current_state%vis_coefficient%data, dim, target_index, current_page(neighbour_location))
537 
538  current_page(neighbour_location)=current_page(neighbour_location)+1
Here is the caller graph for this function:

◆ copy_halo_buffer_to_vis_corners()

subroutine viscosity_mod::copy_halo_buffer_to_vis_corners ( 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_target_index,
integer, intent(in)  y_target_index,
integer, intent(in)  neighbour_location,
integer, dimension(:), intent(inout)  current_page,
type(field_data_wrapper_type), dimension(:), intent(in), optional  source_data 
)
private

Copies the corner halo buffer to the viscosity coefficient field corners.

Parameters
current_stateThe current model state
neighbour_descriptionThe halo swapping description of the neighbour we are accessing the buffer of
corner_locThe corner location
x_target_indexThe X target index for the dimension we are receiving for
y_target_indexThe Y target index for the dimension we are receiving for
neighbour_locationThe location in the local neighbour data stores of this neighbour
current_pageThe current, next, halo swap page to read from (all previous have been read and copied already)
source_dataOptional source data which is written into

Definition at line 550 of file viscosity.F90.

552  type(model_state_type), intent(inout) :: current_state
553  integer, intent(in) :: corner_loc, x_target_index, y_target_index, neighbour_location
554  integer, intent(inout) :: current_page(:)
555  type(neighbour_description_type), intent(inout) :: neighbour_description
556  type(field_data_wrapper_type), dimension(:), intent(in), optional :: source_data
557 
558  call copy_buffer_to_corner(current_state%local_grid, neighbour_description%recv_corner_buffer, &
559  current_state%vis_coefficient%data, corner_loc, x_target_index, y_target_index, current_page(neighbour_location))
560 
561  current_page(neighbour_location)=current_page(neighbour_location)+1
Here is the caller graph for this function:

◆ field_information_retrieval_callback()

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

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

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

Definition at line 72 of file viscosity.F90.

73  type(model_state_type), target, intent(inout) :: current_state
74  character(len=*), intent(in) :: name
75  type(component_field_information_type), intent(out) :: field_information
76  integer :: strcomp
77 
78  ! Field description is the same regardless of the specific field being retrieved
79  field_information%field_type=component_array_field_type
80  field_information%data_type=component_double_data_type
81  field_information%number_dimensions=1
82  field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
83  field_information%enabled=.true.
84 
85  ! Field information for 3d
86  strcomp=index(name, "viscosity_3d_local")
87  if (strcomp .ne. 0) then
88  field_information%field_type=component_array_field_type
89  field_information%number_dimensions=3
90  field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
91  field_information%dimension_sizes(2)=current_state%local_grid%size(y_index)
92  field_information%dimension_sizes(3)=current_state%local_grid%size(x_index)
93  field_information%data_type=component_double_data_type
94 
95  if (name .eq. "tend_u_viscosity_3d_local") then
96  field_information%enabled=l_tend_3d_u
97  else if (name .eq. "tend_v_viscosity_3d_local") then
98  field_information%enabled=l_tend_3d_v
99  else if (name .eq. "tend_w_viscosity_3d_local") then
100  field_information%enabled=l_tend_3d_w
101  else
102  field_information%enabled=.true.
103  end if
104 
105  end if !end 3d check
106 
107  ! Field information for profiles
108  strcomp=index(name, "viscosity_profile_total_local")
109  if (strcomp .ne. 0) then
110  field_information%field_type=component_array_field_type
111  field_information%number_dimensions=1
112  field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
113  field_information%data_type=component_double_data_type
114 
115  if (name .eq. "tend_u_viscosity_profile_total_local") then
116  field_information%enabled=l_tend_pr_tot_u
117  else if (name .eq. "tend_v_viscosity_profile_total_local") then
118  field_information%enabled=l_tend_pr_tot_v
119  else if (name .eq. "tend_w_viscosity_profile_total_local") then
120  field_information%enabled=l_tend_pr_tot_w
121  else
122  field_information%enabled=.true.
123  end if
124 
125  end if !end profile check
126 
Here is the caller graph for this function:

◆ field_value_retrieval_callback()

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

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

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

Definition at line 134 of file viscosity.F90.

135  type(model_state_type), target, intent(inout) :: current_state
136  character(len=*), intent(in) :: name
137  type(component_field_value_type), intent(out) :: field_value
138 
139  if (name .eq. "u_viscosity") then
140  allocate(field_value%real_1d_array(size(u_viscosity)), source=u_viscosity)
141  else if (name .eq. "v_viscosity") then
142  allocate(field_value%real_1d_array(size(v_viscosity)), source=v_viscosity)
143  else if (name .eq. "w_viscosity") then
144  allocate(field_value%real_1d_array(size(w_viscosity)), source=w_viscosity)
145 
146  ! 3d Tendency Fields
147  else if (name .eq. "tend_u_viscosity_3d_local" .and. allocated(tend_3d_u)) then
148  call set_published_field_value(field_value, real_3d_field=tend_3d_u)
149  else if (name .eq. "tend_v_viscosity_3d_local" .and. allocated(tend_3d_v)) then
150  call set_published_field_value(field_value, real_3d_field=tend_3d_v)
151  else if (name .eq. "tend_w_viscosity_3d_local" .and. allocated(tend_3d_w)) then
152  call set_published_field_value(field_value, real_3d_field=tend_3d_w)
153 
154  ! Profile Tendency Fields
155  else if (name .eq. "tend_u_viscosity_profile_total_local" .and. allocated(tend_pr_tot_u)) then
156  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_u)
157  else if (name .eq. "tend_v_viscosity_profile_total_local" .and. allocated(tend_pr_tot_v)) then
158  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_v)
159  else if (name .eq. "tend_w_viscosity_profile_total_local" .and. allocated(tend_pr_tot_w)) then
160  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_w)
161  end if
162 
163 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ finalisation_callback()

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

Definition at line 223 of file viscosity.F90.

224  type(model_state_type), target, intent(inout) :: current_state
225 
226  if (allocated(u_viscosity)) deallocate(u_viscosity)
227  if (allocated(v_viscosity)) deallocate(v_viscosity)
228  if (allocated(w_viscosity)) deallocate(w_viscosity)
229 
230  if (allocated(tend_3d_u)) deallocate(tend_3d_u)
231  if (allocated(tend_3d_v)) deallocate(tend_3d_v)
232  if (allocated(tend_3d_w)) deallocate(tend_3d_w)
233 
234  if (allocated(tend_pr_tot_u)) deallocate(tend_pr_tot_u)
235  if (allocated(tend_pr_tot_v)) deallocate(tend_pr_tot_v)
236  if (allocated(tend_pr_tot_w)) deallocate(tend_pr_tot_w)
237 
Here is the caller graph for this function:

◆ initialisation_callback()

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

Sets up the stencil_mod (used in interpolation) and allocates data for the flux fields.

Parameters
current_stateThe current model state_mod

Definition at line 168 of file viscosity.F90.

169  type(model_state_type), target, intent(inout) :: current_state
170 
171  integer :: z_size, y_size, x_size
172 
173  z_size=current_state%local_grid%size(z_index) + current_state%local_grid%halo_size(z_index) * 2
174  y_size=current_state%local_grid%size(y_index) + current_state%local_grid%halo_size(y_index) * 2
175  x_size=current_state%local_grid%size(x_index) + current_state%local_grid%halo_size(x_index) * 2
176  allocate(current_state%vis_coefficient%data(z_size, y_size, x_size))
177 
178  z_size=current_state%global_grid%size(z_index)
179  allocate(u_viscosity(z_size), v_viscosity(z_size), w_viscosity(z_size))
180 
181  ! Tendency Logicals
182  l_tend_pr_tot_u = current_state%u%active
183  l_tend_pr_tot_v = current_state%v%active
184  l_tend_pr_tot_w = current_state%w%active
185 
186  l_tend_3d_u = current_state%u%active .or. l_tend_pr_tot_u
187  l_tend_3d_v = current_state%v%active .or. l_tend_pr_tot_v
188  l_tend_3d_w = current_state%w%active .or. l_tend_pr_tot_w
189 
190  ! Allocate 3d tendency fields upon availability
191  if (l_tend_3d_u) then
192  allocate( tend_3d_u(current_state%local_grid%size(z_index), &
193  current_state%local_grid%size(y_index), &
194  current_state%local_grid%size(x_index) ) )
195  endif
196  if (l_tend_3d_v) then
197  allocate( tend_3d_v(current_state%local_grid%size(z_index), &
198  current_state%local_grid%size(y_index), &
199  current_state%local_grid%size(x_index) ) )
200  endif
201  if (l_tend_3d_w) then
202  allocate( tend_3d_w(current_state%local_grid%size(z_index), &
203  current_state%local_grid%size(y_index), &
204  current_state%local_grid%size(x_index) ) )
205  endif
206 
207  ! Allocate profile tendency fields upon availability
208  if (l_tend_pr_tot_u) then
209  allocate( tend_pr_tot_u(current_state%local_grid%size(z_index)) )
210  endif
211  if (l_tend_pr_tot_v) then
212  allocate( tend_pr_tot_v(current_state%local_grid%size(z_index)) )
213  endif
214  if (l_tend_pr_tot_w) then
215  allocate( tend_pr_tot_w(current_state%local_grid%size(z_index)) )
216  endif
217 
218  ! Save the sampling_frequency to force diagnostic calculation on select time steps
219  diagnostic_generation_frequency=options_get_integer(current_state%options_database, "sampling_frequency")
220 
Here is the caller graph for this function:

◆ perform_local_data_copy_for_vis()

subroutine viscosity_mod::perform_local_data_copy_for_vis ( type(model_state_type), intent(inout)  current_state,
integer, intent(in)  halo_depth,
logical, intent(in)  involve_corners,
type(field_data_wrapper_type), dimension(:), intent(in), optional  source_data 
)
private

Does local data copying for viscosity coefficient variable halo swap.

Parameters
current_stateThe current model state_mod
source_dataOptional source data which is written into

Definition at line 509 of file viscosity.F90.

510  type(model_state_type), intent(inout) :: current_state
511  integer, intent(in) :: halo_depth
512  logical, intent(in) :: involve_corners
513  type(field_data_wrapper_type), dimension(:), intent(in), optional :: source_data
514 
515  call perform_local_data_copy_for_field(current_state%vis_coefficient%data, current_state%local_grid, &
516  current_state%parallel%my_rank, halo_depth, involve_corners)
Here is the caller graph for this function:

◆ save_precomponent_tendencies()

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

Save the 3d tendencies coming into the component.

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

Definition at line 570 of file viscosity.F90.

571  type(model_state_type), target, intent(in) :: current_state
572  integer, intent(in) :: cxn, cyn, txn, tyn
573 
574  ! Save 3d tendency fields upon request (of 3d or profiles) and availability
575  if (l_tend_3d_u) then
576  tend_3d_u(:,tyn,txn)=current_state%su%data(:,cyn,cxn)
577  endif
578  if (l_tend_3d_v) then
579  tend_3d_v(:,tyn,txn)=current_state%sv%data(:,cyn,cxn)
580  endif
581  if (l_tend_3d_w) then
582  tend_3d_w(:,tyn,txn)=current_state%sw%data(:,cyn,cxn)
583  endif
584 
Here is the caller graph for this function:

◆ set_published_field_value()

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

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

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

Definition at line 627 of file viscosity.F90.

628  type(component_field_value_type), intent(inout) :: field_value
629  real(kind=default_precision), dimension(:), optional :: real_1d_field
630  real(kind=default_precision), dimension(:,:), optional :: real_2d_field
631  real(kind=default_precision), dimension(:,:,:), optional :: real_3d_field
632 
633  if (present(real_1d_field)) then
634  allocate(field_value%real_1d_array(size(real_1d_field)), source=real_1d_field)
635  else if (present(real_2d_field)) then
636  allocate(field_value%real_2d_array(size(real_2d_field, 1), size(real_2d_field, 2)), source=real_2d_field)
637  else if (present(real_3d_field)) then
638  allocate(field_value%real_3d_array(size(real_3d_field, 1), size(real_3d_field, 2), size(real_3d_field, 3)), &
639  source=real_3d_field)
640  end if
Here is the caller graph for this function:

◆ timestep_callback()

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

At each timestep will compute the viscosity U,V,W source terms.

Parameters
current_stateThe current model state

Definition at line 242 of file viscosity.F90.

243  type(model_state_type), target, intent(inout) :: current_state
244 
245  integer :: local_y, locaL_x, k, target_x_index, target_y_index
246  real(kind=default_precision), dimension(current_state%local_grid%size(Z_INDEX)) :: tau12, tau12_ym1, tau12m1, &
247  tau11, tau22, tau22_yp1, tau33, tau23_ym1, tau11p1, tau13, tau13m1, tau23
248 
249  local_y=current_state%column_local_y
250  local_x=current_state%column_local_x
251  target_y_index=local_y-current_state%local_grid%halo_size(y_index)
252  target_x_index=local_x-current_state%local_grid%halo_size(x_index)
253 
254  ! Zero profile tendency totals on first instance in the sum
255  if (current_state%first_timestep_column) then
256  if (l_tend_pr_tot_u) then
257  tend_pr_tot_u(:)= 0.0_default_precision
258  endif
259  if (l_tend_pr_tot_v) then
260  tend_pr_tot_v(:)= 0.0_default_precision
261  endif
262  if (l_tend_pr_tot_w) then
263  tend_pr_tot_w(:)= 0.0_default_precision
264  endif
265  endif ! zero totals
266 
267  if (.not. current_state%use_viscosity_and_diffusion .or. current_state%halo_column) return
268 
269  if (current_state%viscosity_halo_swap_state%swap_in_progress) then
270  ! If there is a viscosity halo swap in progress then complete it
271  call complete_nonblocking_halo_swap(current_state, current_state%viscosity_halo_swap_state, &
272  perform_local_data_copy_for_vis, copy_halo_buffer_to_vis, copy_halo_buffer_to_vis_corners)
273  end if
274 
275  if (mod(current_state%timestep, diagnostic_generation_frequency) == 0) then
276  call save_precomponent_tendencies(current_state, local_x, local_y, target_x_index, target_y_index)
277  end if
278 
279  if (current_state%field_stepping == forward_stepping) then
280  call calculate_tau(current_state, local_y, local_x, current_state%u, current_state%v, current_state%w, tau12, tau12_ym1, &
281  tau12m1, tau11, tau22, tau22_yp1, tau33, tau11p1, tau13, tau13m1, tau23, tau23_ym1)
282  else
283  call calculate_tau(current_state, local_y, local_x, current_state%zu, current_state%zv, current_state%zw, tau12, tau12_ym1, &
284  tau12m1, tau11, tau22, tau22_yp1, tau33, tau11p1, tau13, tau13m1, tau23, tau23_ym1)
285  end if
286  call calculate_viscous_sources(current_state, local_y, local_x, tau12, tau12_ym1, tau12m1, tau11, tau22, tau22_yp1, tau33, &
287  tau11p1, tau13, tau13m1, tau23, tau23_ym1)
288 
289  if (mod(current_state%timestep, diagnostic_generation_frequency) == 0) then
290  call compute_component_tendencies(current_state, local_x, local_y, target_x_index, target_y_index)
291  end if
292 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ viscosity_get_descriptor()

type(component_descriptor_type) function, public viscosity_mod::viscosity_get_descriptor

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

Returns
The termination check component descriptor

Definition at line 43 of file viscosity.F90.

44  viscosity_get_descriptor%name="viscosity"
45  viscosity_get_descriptor%version=0.1
46  viscosity_get_descriptor%initialisation=>initialisation_callback
47  viscosity_get_descriptor%timestep=>timestep_callback
48  viscosity_get_descriptor%finalisation=>finalisation_callback
49 
50  viscosity_get_descriptor%field_value_retrieval=>field_value_retrieval_callback
51  viscosity_get_descriptor%field_information_retrieval=>field_information_retrieval_callback
52  allocate(viscosity_get_descriptor%published_fields(3+3+3))
53 
54  viscosity_get_descriptor%published_fields(1)="u_viscosity"
55  viscosity_get_descriptor%published_fields(2)="v_viscosity"
56  viscosity_get_descriptor%published_fields(3)="w_viscosity"
57 
58  viscosity_get_descriptor%published_fields(3+1)= "tend_u_viscosity_3d_local"
59  viscosity_get_descriptor%published_fields(3+2)= "tend_v_viscosity_3d_local"
60  viscosity_get_descriptor%published_fields(3+3)= "tend_w_viscosity_3d_local"
61 
62  viscosity_get_descriptor%published_fields(3+3+1)= "tend_u_viscosity_profile_total_local"
63  viscosity_get_descriptor%published_fields(3+3+2)= "tend_v_viscosity_profile_total_local"
64  viscosity_get_descriptor%published_fields(3+3+3)= "tend_w_viscosity_profile_total_local"
65 
Here is the call graph for this function:

Variable Documentation

◆ diagnostic_generation_frequency

integer viscosity_mod::diagnostic_generation_frequency
private

Definition at line 35 of file viscosity.F90.

35  integer :: diagnostic_generation_frequency

◆ l_tend_3d_u

logical viscosity_mod::l_tend_3d_u
private

Definition at line 29 of file viscosity.F90.

29  logical :: l_tend_3d_u, l_tend_3d_v, l_tend_3d_w

◆ l_tend_3d_v

logical viscosity_mod::l_tend_3d_v
private

Definition at line 29 of file viscosity.F90.

◆ l_tend_3d_w

logical viscosity_mod::l_tend_3d_w
private

Definition at line 29 of file viscosity.F90.

◆ l_tend_pr_tot_u

logical viscosity_mod::l_tend_pr_tot_u
private

Definition at line 33 of file viscosity.F90.

33  logical :: l_tend_pr_tot_u, l_tend_pr_tot_v, l_tend_pr_tot_w

◆ l_tend_pr_tot_v

logical viscosity_mod::l_tend_pr_tot_v
private

Definition at line 33 of file viscosity.F90.

◆ l_tend_pr_tot_w

logical viscosity_mod::l_tend_pr_tot_w
private

Definition at line 33 of file viscosity.F90.

◆ tend_3d_u

real(kind=default_precision), dimension(:,:,:), allocatable viscosity_mod::tend_3d_u
private

Definition at line 27 of file viscosity.F90.

27  real(kind=default_precision), dimension(:,:,:), allocatable :: &
28  tend_3d_u, tend_3d_v, tend_3d_w

◆ tend_3d_v

real(kind=default_precision), dimension(:,:,:), allocatable viscosity_mod::tend_3d_v
private

Definition at line 27 of file viscosity.F90.

◆ tend_3d_w

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

Definition at line 27 of file viscosity.F90.

◆ tend_pr_tot_u

real(kind=default_precision), dimension(:), allocatable viscosity_mod::tend_pr_tot_u
private

Definition at line 31 of file viscosity.F90.

31  real(kind=default_precision), dimension(:), allocatable :: &
32  tend_pr_tot_u, tend_pr_tot_v, tend_pr_tot_w

◆ tend_pr_tot_v

real(kind=default_precision), dimension(:), allocatable viscosity_mod::tend_pr_tot_v
private

Definition at line 31 of file viscosity.F90.

◆ tend_pr_tot_w

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

Definition at line 31 of file viscosity.F90.

◆ u_viscosity

real(kind=default_precision), dimension(:), allocatable viscosity_mod::u_viscosity
private

Definition at line 23 of file viscosity.F90.

23  real(kind=default_precision), dimension(:), allocatable :: u_viscosity, v_viscosity, w_viscosity

◆ v_viscosity

real(kind=default_precision), dimension(:), allocatable viscosity_mod::v_viscosity
private

Definition at line 23 of file viscosity.F90.

◆ w_viscosity

real(kind=default_precision), dimension(:), allocatable viscosity_mod::w_viscosity
private

Definition at line 23 of file viscosity.F90.

grids_mod::x_index
integer, parameter, public x_index
Definition: grids.F90:14
optionsdatabase_mod::options_get_integer
integer function, public options_get_integer(options_database, key, index)
Retrieves an integer value from the database that matches the provided key.
Definition: optionsdatabase.F90:217
grids_mod::y_index
integer, parameter, public y_index
Definition: grids.F90:14
grids_mod::z_index
integer, parameter, public z_index
Grid index parameters.
Definition: grids.F90:14
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