MONC
Functions/Subroutines | Variables
diffusion_mod Module Reference

Diffusion on the TH and Q fields. More...

Functions/Subroutines

type(component_descriptor_type) function, public diffusion_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 diffusion source terms for TH and Q fields per column if these fields are active. More...
 
subroutine perform_q_diffusion (current_state, local_y, local_x)
 Computes the diffusion source terms for each Q field. More...
 
subroutine perform_th_diffusion (current_state, local_y, local_x)
 Computes the diffusion source terms for the theta field. More...
 
subroutine general_diffusion (current_state, field, source_field, local_y, local_x, diagnostics)
 General diffusion computation for any field which is provided as arguments. Works in a column. More...
 
subroutine perform_local_data_copy_for_diff (current_state, halo_depth, involve_corners, source_data)
 Does local data copying for diffusion coefficient variable halo swap. More...
 
subroutine copy_halo_buffer_to_diff (current_state, neighbour_description, dim, target_index, neighbour_location, current_page, source_data)
 Copies the halo buffer to halo location for the diffusion coefficient field. More...
 
subroutine copy_halo_buffer_to_diff_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 diffusion 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 th_diffusion
 
real(kind=default_precision), dimension(:,:), allocatable q_diffusion
 
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_th
 
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_qv
 
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_ql
 
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_qi
 
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_qr
 
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_qs
 
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_qg
 
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_tabs
 
logical l_tend_3d_th
 
logical l_tend_3d_qv
 
logical l_tend_3d_ql
 
logical l_tend_3d_qi
 
logical l_tend_3d_qr
 
logical l_tend_3d_qs
 
logical l_tend_3d_qg
 
logical l_tend_3d_tabs
 
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_th
 
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_qv
 
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_ql
 
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_qi
 
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_qr
 
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_qs
 
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_qg
 
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_tabs
 
logical l_tend_pr_tot_th
 
logical l_tend_pr_tot_qv
 
logical l_tend_pr_tot_ql
 
logical l_tend_pr_tot_qi
 
logical l_tend_pr_tot_qr
 
logical l_tend_pr_tot_qs
 
logical l_tend_pr_tot_qg
 
logical l_tend_pr_tot_tabs
 
integer iqv =0
 
integer iql =0
 
integer iqr =0
 
integer iqi =0
 
integer iqs =0
 
integer iqg =0
 
integer diagnostic_generation_frequency
 

Detailed Description

Diffusion on the TH and Q fields.

Function/Subroutine Documentation

◆ compute_component_tendencies()

subroutine diffusion_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 626 of file diffusion.F90.

627  type(model_state_type), target, intent(inout) :: current_state
628  integer, intent(in) :: cxn, cyn, txn, tyn
629 
630  ! Calculate change in tendency due to component
631  if (l_tend_3d_th) then
632  tend_3d_th(:,tyn,txn)=current_state%sth%data(:,cyn,cxn) - tend_3d_th(:,tyn,txn)
633  endif
634  if (l_tend_3d_qv) then
635  tend_3d_qv(:,tyn,txn)=current_state%sq(iqv)%data(:,cyn,cxn) - tend_3d_qv(:,tyn,txn)
636  endif
637  if (l_tend_3d_ql) then
638  tend_3d_ql(:,tyn,txn)=current_state%sq(iql)%data(:,cyn,cxn) - tend_3d_ql(:,tyn,txn)
639  endif
640  if (l_tend_3d_qi) then
641  tend_3d_qi(:,tyn,txn)=current_state%sq(iqi)%data(:,cyn,cxn) - tend_3d_qi(:,tyn,txn)
642  endif
643  if (l_tend_3d_qr) then
644  tend_3d_qr(:,tyn,txn)=current_state%sq(iqr)%data(:,cyn,cxn) - tend_3d_qr(:,tyn,txn)
645  endif
646  if (l_tend_3d_qs) then
647  tend_3d_qs(:,tyn,txn)=current_state%sq(iqs)%data(:,cyn,cxn) - tend_3d_qs(:,tyn,txn)
648  endif
649  if (l_tend_3d_qg) then
650  tend_3d_qg(:,tyn,txn)=current_state%sq(iqg)%data(:,cyn,cxn) - tend_3d_qg(:,tyn,txn)
651  endif
652  if (l_tend_3d_tabs) then
653  tend_3d_tabs(:,tyn,txn)= &
654  current_state%sth%data(:,cyn,cxn) * current_state%global_grid%configuration%vertical%rprefrcp(:) &
655  - tend_3d_tabs(:,tyn,txn)
656  endif
657 
658  ! Add local tendency fields to the profile total
659  if (l_tend_pr_tot_th) then
660  tend_pr_tot_th(:)=tend_pr_tot_th(:) + tend_3d_th(:,tyn,txn)
661  endif
662  if (l_tend_pr_tot_qv) then
663  tend_pr_tot_qv(:)=tend_pr_tot_qv(:) + tend_3d_qv(:,tyn,txn)
664  endif
665  if (l_tend_pr_tot_ql) then
666  tend_pr_tot_ql(:)=tend_pr_tot_ql(:) + tend_3d_ql(:,tyn,txn)
667  endif
668  if (l_tend_pr_tot_qi) then
669  tend_pr_tot_qi(:)=tend_pr_tot_qi(:) + tend_3d_qi(:,tyn,txn)
670  endif
671  if (l_tend_pr_tot_qr) then
672  tend_pr_tot_qr(:)=tend_pr_tot_qr(:) + tend_3d_qr(:,tyn,txn)
673  endif
674  if (l_tend_pr_tot_qs) then
675  tend_pr_tot_qs(:)=tend_pr_tot_qs(:) + tend_3d_qs(:,tyn,txn)
676  endif
677  if (l_tend_pr_tot_qg) then
678  tend_pr_tot_qg(:)=tend_pr_tot_qg(:) + tend_3d_qg(:,tyn,txn)
679  endif
680  if (l_tend_pr_tot_tabs) then
681  tend_pr_tot_tabs(:)=tend_pr_tot_tabs(:) + tend_3d_tabs(:,tyn,txn)
682  endif
683 
Here is the caller graph for this function:

◆ copy_halo_buffer_to_diff()

subroutine diffusion_mod::copy_halo_buffer_to_diff ( 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 diffusion 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 544 of file diffusion.F90.

546  type(model_state_type), intent(inout) :: current_state
547  integer, intent(in) :: dim, target_index, neighbour_location
548  integer, intent(inout) :: current_page(:)
549  type(neighbour_description_type), intent(inout) :: neighbour_description
550  type(field_data_wrapper_type), dimension(:), intent(in), optional :: source_data
551 
552  call copy_buffer_to_field(current_state%local_grid, neighbour_description%recv_halo_buffer, &
553  current_state%diff_coefficient%data, dim, target_index, current_page(neighbour_location))
554 
555  current_page(neighbour_location)=current_page(neighbour_location)+1
Here is the caller graph for this function:

◆ copy_halo_buffer_to_diff_corners()

subroutine diffusion_mod::copy_halo_buffer_to_diff_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 diffusion 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 567 of file diffusion.F90.

569  type(model_state_type), intent(inout) :: current_state
570  integer, intent(in) :: corner_loc, x_target_index, y_target_index, neighbour_location
571  integer, intent(inout) :: current_page(:)
572  type(neighbour_description_type), intent(inout) :: neighbour_description
573  type(field_data_wrapper_type), dimension(:), intent(in), optional :: source_data
574 
575  call copy_buffer_to_corner(current_state%local_grid, neighbour_description%recv_corner_buffer, &
576  current_state%diff_coefficient%data, corner_loc, x_target_index, y_target_index, current_page(neighbour_location))
577 
578  current_page(neighbour_location)=current_page(neighbour_location)+1
Here is the caller graph for this function:

◆ diffusion_get_descriptor()

type(component_descriptor_type) function, public diffusion_mod::diffusion_get_descriptor

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

Returns
The termination check component descriptor

Definition at line 48 of file diffusion.F90.

49  diffusion_get_descriptor%name="diffusion"
50  diffusion_get_descriptor%version=0.1
51  diffusion_get_descriptor%initialisation=>initialisation_callback
52  diffusion_get_descriptor%timestep=>timestep_callback
53  diffusion_get_descriptor%finalisation=>finalisation_callback
54 
55  diffusion_get_descriptor%field_value_retrieval=>field_value_retrieval_callback
56  diffusion_get_descriptor%field_information_retrieval=>field_information_retrieval_callback
57  allocate(diffusion_get_descriptor%published_fields(2+8+8))
58 
59  diffusion_get_descriptor%published_fields(1)="th_diffusion"
60  diffusion_get_descriptor%published_fields(2)="q_diffusion"
61 
62  diffusion_get_descriptor%published_fields(2+1)="tend_th_diffusion_3d_local"
63  diffusion_get_descriptor%published_fields(2+2)="tend_qv_diffusion_3d_local"
64  diffusion_get_descriptor%published_fields(2+3)="tend_ql_diffusion_3d_local"
65  diffusion_get_descriptor%published_fields(2+4)="tend_qi_diffusion_3d_local"
66  diffusion_get_descriptor%published_fields(2+5)="tend_qr_diffusion_3d_local"
67  diffusion_get_descriptor%published_fields(2+6)="tend_qs_diffusion_3d_local"
68  diffusion_get_descriptor%published_fields(2+7)="tend_qg_diffusion_3d_local"
69  diffusion_get_descriptor%published_fields(2+8)="tend_tabs_diffusion_3d_local"
70 
71  diffusion_get_descriptor%published_fields(2+8+1)="tend_th_diffusion_profile_total_local"
72  diffusion_get_descriptor%published_fields(2+8+2)="tend_qv_diffusion_profile_total_local"
73  diffusion_get_descriptor%published_fields(2+8+3)="tend_ql_diffusion_profile_total_local"
74  diffusion_get_descriptor%published_fields(2+8+4)="tend_qi_diffusion_profile_total_local"
75  diffusion_get_descriptor%published_fields(2+8+5)="tend_qr_diffusion_profile_total_local"
76  diffusion_get_descriptor%published_fields(2+8+6)="tend_qs_diffusion_profile_total_local"
77  diffusion_get_descriptor%published_fields(2+8+7)="tend_qg_diffusion_profile_total_local"
78  diffusion_get_descriptor%published_fields(2+8+8)="tend_tabs_diffusion_profile_total_local"
79 
Here is the call graph for this function:

◆ field_information_retrieval_callback()

subroutine diffusion_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 87 of file diffusion.F90.

88  type(model_state_type), target, intent(inout) :: current_state
89  character(len=*), intent(in) :: name
90  type(component_field_information_type), intent(out) :: field_information
91  integer :: strcomp
92 
93  ! Field description is the same regardless of the specific field being retrieved
94  field_information%field_type=component_array_field_type
95  field_information%data_type=component_double_data_type
96  if (name .eq. "q_diffusion") then
97  field_information%number_dimensions=2
98  else
99  field_information%number_dimensions=1
100  end if
101  field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
102  if (name .eq. "q_diffusion") field_information%dimension_sizes(2)=current_state%number_q_fields
103  field_information%enabled=.true.
104 
105  ! Field information for 3d
106  strcomp=index(name, "diffusion_3d_local")
107  if (strcomp .ne. 0) then
108  field_information%field_type=component_array_field_type
109  field_information%number_dimensions=3
110  field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
111  field_information%dimension_sizes(2)=current_state%local_grid%size(y_index)
112  field_information%dimension_sizes(3)=current_state%local_grid%size(x_index)
113  field_information%data_type=component_double_data_type
114 
115  if (name .eq. "tend_th_diffusion_3d_local") then
116  field_information%enabled=l_tend_3d_th
117  else if (name .eq. "tend_qv_diffusion_3d_local") then
118  field_information%enabled=l_tend_3d_qv
119  else if (name .eq. "tend_ql_diffusion_3d_local") then
120  field_information%enabled=l_tend_3d_ql
121  else if (name .eq. "tend_qi_diffusion_3d_local") then
122  field_information%enabled=l_tend_3d_qi
123  else if (name .eq. "tend_qr_diffusion_3d_local") then
124  field_information%enabled=l_tend_3d_qr
125  else if (name .eq. "tend_qs_diffusion_3d_local") then
126  field_information%enabled=l_tend_3d_qs
127  else if (name .eq. "tend_qg_diffusion_3d_local") then
128  field_information%enabled=l_tend_3d_qg
129  else if (name .eq. "tend_tabs_diffusion_3d_local") then
130  field_information%enabled=l_tend_3d_tabs
131  else
132  field_information%enabled=.true.
133  end if
134 
135  end if !end 3d check
136 
137  ! Field information for profiles
138  strcomp=index(name, "diffusion_profile_total_local")
139  if (strcomp .ne. 0) then
140  field_information%field_type=component_array_field_type
141  field_information%number_dimensions=1
142  field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
143  field_information%data_type=component_double_data_type
144 
145  if (name .eq. "tend_th_diffusion_profile_total_local") then
146  field_information%enabled=l_tend_pr_tot_th
147  else if (name .eq. "tend_qv_diffusion_profile_total_local") then
148  field_information%enabled=l_tend_pr_tot_qv
149  else if (name .eq. "tend_ql_diffusion_profile_total_local") then
150  field_information%enabled=l_tend_pr_tot_ql
151  else if (name .eq. "tend_qi_diffusion_profile_total_local") then
152  field_information%enabled=l_tend_pr_tot_qi
153  else if (name .eq. "tend_qr_diffusion_profile_total_local") then
154  field_information%enabled=l_tend_pr_tot_qr
155  else if (name .eq. "tend_qs_diffusion_profile_total_local") then
156  field_information%enabled=l_tend_pr_tot_qs
157  else if (name .eq. "tend_qg_diffusion_profile_total_local") then
158  field_information%enabled=l_tend_pr_tot_qg
159  else if (name .eq. "tend_tabs_diffusion_profile_total_local") then
160  field_information%enabled=l_tend_pr_tot_tabs
161  else
162  field_information%enabled=.true.
163  end if
164 
165  end if !end profile check
166 
Here is the caller graph for this function:

◆ field_value_retrieval_callback()

subroutine diffusion_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 173 of file diffusion.F90.

174  type(model_state_type), target, intent(inout) :: current_state
175  character(len=*), intent(in) :: name
176  type(component_field_value_type), intent(out) :: field_value
177 
178  if (name .eq. "th_diffusion") then
179  call set_published_field_value(field_value, real_1d_field=th_diffusion)
180  else if (name .eq. "q_diffusion") then
181  call set_published_field_value(field_value, real_2d_field=q_diffusion)
182 
183  ! 3d Tendency Fields
184  else if (name .eq. "tend_th_diffusion_3d_local" .and. allocated(tend_3d_th)) then
185  call set_published_field_value(field_value, real_3d_field=tend_3d_th)
186  else if (name .eq. "tend_qv_diffusion_3d_local" .and. allocated(tend_3d_qv)) then
187  call set_published_field_value(field_value, real_3d_field=tend_3d_qv)
188  else if (name .eq. "tend_ql_diffusion_3d_local" .and. allocated(tend_3d_ql)) then
189  call set_published_field_value(field_value, real_3d_field=tend_3d_ql)
190  else if (name .eq. "tend_qi_diffusion_3d_local" .and. allocated(tend_3d_qi)) then
191  call set_published_field_value(field_value, real_3d_field=tend_3d_qi)
192  else if (name .eq. "tend_qr_diffusion_3d_local" .and. allocated(tend_3d_qr)) then
193  call set_published_field_value(field_value, real_3d_field=tend_3d_qr)
194  else if (name .eq. "tend_qs_diffusion_3d_local" .and. allocated(tend_3d_qs)) then
195  call set_published_field_value(field_value, real_3d_field=tend_3d_qs)
196  else if (name .eq. "tend_qg_diffusion_3d_local" .and. allocated(tend_3d_qg)) then
197  call set_published_field_value(field_value, real_3d_field=tend_3d_qg)
198  else if (name .eq. "tend_tabs_diffusion_3d_local" .and. allocated(tend_3d_tabs)) then
199  call set_published_field_value(field_value, real_3d_field=tend_3d_tabs)
200 
201  ! Profile Tendency Fields
202  else if (name .eq. "tend_th_diffusion_profile_total_local" .and. allocated(tend_pr_tot_th)) then
203  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_th)
204  else if (name .eq. "tend_qv_diffusion_profile_total_local" .and. allocated(tend_pr_tot_qv)) then
205  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_qv)
206  else if (name .eq. "tend_ql_diffusion_profile_total_local" .and. allocated(tend_pr_tot_ql)) then
207  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_ql)
208  else if (name .eq. "tend_qi_diffusion_profile_total_local" .and. allocated(tend_pr_tot_qi)) then
209  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_qi)
210  else if (name .eq. "tend_qr_diffusion_profile_total_local" .and. allocated(tend_pr_tot_qr)) then
211  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_qr)
212  else if (name .eq. "tend_qs_diffusion_profile_total_local" .and. allocated(tend_pr_tot_qs)) then
213  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_qs)
214  else if (name .eq. "tend_qg_diffusion_profile_total_local" .and. allocated(tend_pr_tot_qg)) then
215  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_qg)
216  else if (name .eq. "tend_tabs_diffusion_profile_total_local" .and. allocated(tend_pr_tot_tabs)) then
217  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_tabs)
218  end if
219 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ finalisation_callback()

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

Definition at line 342 of file diffusion.F90.

343  type(model_state_type), target, intent(inout) :: current_state
344 
345  if (allocated(th_diffusion)) deallocate(th_diffusion)
346  if (allocated(q_diffusion)) deallocate(q_diffusion)
347 
348  if (allocated(tend_3d_th)) deallocate(tend_3d_th)
349  if (allocated(tend_3d_qv)) deallocate(tend_3d_qv)
350  if (allocated(tend_3d_ql)) deallocate(tend_3d_ql)
351  if (allocated(tend_3d_qi)) deallocate(tend_3d_qi)
352  if (allocated(tend_3d_qr)) deallocate(tend_3d_qr)
353  if (allocated(tend_3d_qs)) deallocate(tend_3d_qs)
354  if (allocated(tend_3d_qg)) deallocate(tend_3d_qg)
355  if (allocated(tend_3d_tabs)) deallocate(tend_3d_tabs)
356 
357  if (allocated(tend_pr_tot_th)) deallocate(tend_pr_tot_th)
358  if (allocated(tend_pr_tot_qv)) deallocate(tend_pr_tot_qv)
359  if (allocated(tend_pr_tot_ql)) deallocate(tend_pr_tot_ql)
360  if (allocated(tend_pr_tot_qi)) deallocate(tend_pr_tot_qi)
361  if (allocated(tend_pr_tot_qr)) deallocate(tend_pr_tot_qr)
362  if (allocated(tend_pr_tot_qs)) deallocate(tend_pr_tot_qs)
363  if (allocated(tend_pr_tot_qg)) deallocate(tend_pr_tot_qg)
364  if (allocated(tend_pr_tot_tabs)) deallocate(tend_pr_tot_tabs)
365 
Here is the caller graph for this function:

◆ general_diffusion()

subroutine diffusion_mod::general_diffusion ( type(model_state_type), intent(inout), target  current_state,
type(prognostic_field_type), intent(inout)  field,
type(prognostic_field_type), intent(inout)  source_field,
integer, intent(in)  local_y,
integer, intent(in)  local_x,
real(kind=default_precision), dimension(:), intent(inout), optional  diagnostics 
)
private

General diffusion computation for any field which is provided as arguments. Works in a column.

Parameters
current_stateThe current model state
fieldThe field to take values from, typically zth or zq(n)
source_fieldThe source target field to update, typically sth or sq(n)
local_yLocal Y index
local_xLocal X index

Definition at line 477 of file diffusion.F90.

478  type(model_state_type), target, intent(inout) :: current_state
479  type(prognostic_field_type), intent(inout) :: field, source_field
480  integer, intent(in) :: local_y, local_x
481  real(kind=default_precision), dimension(:), intent(inout), optional :: diagnostics
482 
483  real(kind=default_precision) :: term
484  integer :: k
485 
486  do k=2, current_state%local_grid%size(z_index)
487  term=current_state%global_grid%configuration%horizontal%cx2*0.25_default_precision*&
488  (((current_state%diff_coefficient%data(k, local_y, local_x)+&
489  current_state%diff_coefficient%data(k, local_y, local_x-1))+&
490  (current_state%diff_coefficient%data(k-1, local_y, local_x)+&
491  current_state%diff_coefficient%data(k-1, local_y, local_x-1)))&
492  *(field%data(k, local_y, local_x-1)-field%data(k, local_y, local_x)) -&
493  ((current_state%diff_coefficient%data(k, local_y, local_x+1)+&
494  current_state%diff_coefficient%data(k, local_y, local_x))+&
495  (current_state%diff_coefficient%data(k-1, local_y, local_x+1)+&
496  current_state%diff_coefficient%data(k-1, local_y, local_x)))&
497  *(field%data(k, local_y, local_x)-field%data(k, local_y, local_x+1)) )&
498  +current_state%global_grid%configuration%horizontal%cy2*0.25_default_precision*(&
499  ((current_state%diff_coefficient%data(k, local_y, local_x)+&
500  current_state%diff_coefficient%data(k, local_y-1, local_x))+&
501  (current_state%diff_coefficient%data(k-1, local_y, local_x)+&
502  current_state%diff_coefficient%data(k-1, local_y-1, local_x)))&
503  *(field%data(k, local_y-1, local_x)-field%data(k, local_y, local_x)) -&
504  ((current_state%diff_coefficient%data(k, local_y+1, local_x)+&
505  current_state%diff_coefficient%data(k, local_y, local_x))+&
506  (current_state%diff_coefficient%data(k-1, local_y+1, local_x)+&
507  current_state%diff_coefficient%data(k-1, local_y, local_x)))&
508  *(field%data(k, local_y, local_x)-field%data(k, local_y+1, local_x)) )&
509  +( current_state%global_grid%configuration%vertical%czb(k)*&
510  current_state%diff_coefficient%data(k-1, local_y, local_x)*&
511  (field%data(k-1, local_y, local_x)-field%data(k, local_y, local_x)))
512 
513  if (k .lt. current_state%local_grid%size(z_index)) then
514  term=term - current_state%global_grid%configuration%vertical%cza(k)*&
515  current_state%diff_coefficient%data(k, local_y, local_x)*&
516  (field%data(k, local_y, local_x)-field%data(k+1, local_y, local_x))
517  end if
518  source_field%data(k, local_y, local_x)=source_field%data(k, local_y, local_x)+term
519  if (present(diagnostics)) diagnostics(k)=term
520  end do
Here is the caller graph for this function:

◆ initialisation_callback()

subroutine diffusion_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 224 of file diffusion.F90.

225  type(model_state_type), target, intent(inout) :: current_state
226 
227  integer :: z_size, y_size, x_size
228  logical :: l_qdiag
229 
230  z_size=current_state%local_grid%size(z_index) + current_state%local_grid%halo_size(z_index) * 2
231  y_size=current_state%local_grid%size(y_index) + current_state%local_grid%halo_size(y_index) * 2
232  x_size=current_state%local_grid%size(x_index) + current_state%local_grid%halo_size(x_index) * 2
233  allocate(current_state%diff_coefficient%data(z_size, y_size, x_size))
234 
235  z_size=current_state%global_grid%size(z_index)
236  allocate(th_diffusion(z_size))
237  allocate(q_diffusion(z_size, current_state%number_q_fields))
238 
239  ! Set tendency diagnostic logicals based on availability
240  ! Need to use 3d tendencies to compute the profiles, so they will be allocated
241  ! in the case where profiles are available
242  l_qdiag = (.not. current_state%passive_q .and. current_state%number_q_fields .gt. 0)
243 
244  l_tend_pr_tot_th = current_state%th%active
245  l_tend_pr_tot_qv = l_qdiag .and. current_state%number_q_fields .ge. 1
246  l_tend_pr_tot_ql = l_qdiag .and. current_state%number_q_fields .ge. 2
247  l_tend_pr_tot_qi = l_qdiag .and. current_state%number_q_fields .ge. 11
248  l_tend_pr_tot_qr = l_qdiag .and. current_state%number_q_fields .ge. 11
249  l_tend_pr_tot_qs = l_qdiag .and. current_state%number_q_fields .ge. 11
250  l_tend_pr_tot_qg = l_qdiag .and. current_state%number_q_fields .ge. 11
251  l_tend_pr_tot_tabs = l_tend_pr_tot_th
252 
253  l_tend_3d_th = current_state%th%active .or. l_tend_pr_tot_th
254  l_tend_3d_qv = (l_qdiag .and. current_state%number_q_fields .ge. 1) .or. l_tend_pr_tot_qv
255  l_tend_3d_ql = (l_qdiag .and. current_state%number_q_fields .ge. 2) .or. l_tend_pr_tot_ql
256  l_tend_3d_qi = (l_qdiag .and. current_state%number_q_fields .ge. 11) .or. l_tend_pr_tot_qi
257  l_tend_3d_qr = (l_qdiag .and. current_state%number_q_fields .ge. 11) .or. l_tend_pr_tot_qr
258  l_tend_3d_qs = (l_qdiag .and. current_state%number_q_fields .ge. 11) .or. l_tend_pr_tot_qs
259  l_tend_3d_qg = (l_qdiag .and. current_state%number_q_fields .ge. 11) .or. l_tend_pr_tot_qg
260  l_tend_3d_tabs = l_tend_3d_th
261 
262  ! Allocate 3d tendency fields upon availability
263  if (l_tend_3d_th) then
264  allocate( tend_3d_th(current_state%local_grid%size(z_index), &
265  current_state%local_grid%size(y_index), &
266  current_state%local_grid%size(x_index) ) )
267  endif
268  if (l_tend_3d_qv) then
269  iqv=get_q_index(standard_q_names%VAPOUR, 'diffusion')
270  allocate( tend_3d_qv(current_state%local_grid%size(z_index), &
271  current_state%local_grid%size(y_index), &
272  current_state%local_grid%size(x_index) ) )
273  endif
274  if (l_tend_3d_ql) then
275  iql=get_q_index(standard_q_names%CLOUD_LIQUID_MASS, 'diffusion')
276  allocate( tend_3d_ql(current_state%local_grid%size(z_index), &
277  current_state%local_grid%size(y_index), &
278  current_state%local_grid%size(x_index) ) )
279  endif
280  if (l_tend_3d_qi) then
281  iqi=get_q_index(standard_q_names%ICE_MASS, 'diffusion')
282  allocate( tend_3d_qi(current_state%local_grid%size(z_index), &
283  current_state%local_grid%size(y_index), &
284  current_state%local_grid%size(x_index) ) )
285  endif
286  if (l_tend_3d_qr) then
287  iqr=get_q_index(standard_q_names%RAIN_MASS, 'diffusion')
288  allocate( tend_3d_qr(current_state%local_grid%size(z_index), &
289  current_state%local_grid%size(y_index), &
290  current_state%local_grid%size(x_index) ) )
291  endif
292  if (l_tend_3d_qs) then
293  iqs=get_q_index(standard_q_names%SNOW_MASS, 'diffusion')
294  allocate( tend_3d_qs(current_state%local_grid%size(z_index), &
295  current_state%local_grid%size(y_index), &
296  current_state%local_grid%size(x_index) ) )
297  endif
298  if (l_tend_3d_qg) then
299  iqg=get_q_index(standard_q_names%GRAUPEL_MASS, 'diffusion')
300  allocate( tend_3d_qg(current_state%local_grid%size(z_index), &
301  current_state%local_grid%size(y_index), &
302  current_state%local_grid%size(x_index) ) )
303  endif
304  if (l_tend_3d_tabs) then
305  allocate( tend_3d_tabs(current_state%local_grid%size(z_index), &
306  current_state%local_grid%size(y_index), &
307  current_state%local_grid%size(x_index) ) )
308  endif
309 
310  ! Allocate profile tendency fields upon availability
311  if (l_tend_pr_tot_th) then
312  allocate( tend_pr_tot_th(current_state%local_grid%size(z_index)) )
313  endif
314  if (l_tend_pr_tot_qv) then
315  allocate( tend_pr_tot_qv(current_state%local_grid%size(z_index)) )
316  endif
317  if (l_tend_pr_tot_ql) then
318  allocate( tend_pr_tot_ql(current_state%local_grid%size(z_index)) )
319  endif
320  if (l_tend_pr_tot_qi) then
321  allocate( tend_pr_tot_qi(current_state%local_grid%size(z_index)) )
322  endif
323  if (l_tend_pr_tot_qr) then
324  allocate( tend_pr_tot_qr(current_state%local_grid%size(z_index)) )
325  endif
326  if (l_tend_pr_tot_qs) then
327  allocate( tend_pr_tot_qs(current_state%local_grid%size(z_index)) )
328  endif
329  if (l_tend_pr_tot_qg) then
330  allocate( tend_pr_tot_qg(current_state%local_grid%size(z_index)) )
331  endif
332  if (l_tend_pr_tot_tabs) then
333  allocate( tend_pr_tot_tabs(current_state%local_grid%size(z_index)) )
334  endif
335 
336  ! Save the sampling_frequency to force diagnostic calculation on select time steps
337  diagnostic_generation_frequency=options_get_integer(current_state%options_database, "sampling_frequency")
338 
Here is the caller graph for this function:

◆ perform_local_data_copy_for_diff()

subroutine diffusion_mod::perform_local_data_copy_for_diff ( 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 diffusion coefficient variable halo swap.

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

Definition at line 526 of file diffusion.F90.

527  type(model_state_type), intent(inout) :: current_state
528  integer, intent(in) :: halo_depth
529  logical, intent(in) :: involve_corners
530  type(field_data_wrapper_type), dimension(:), intent(in), optional :: source_data
531 
532  call perform_local_data_copy_for_field(current_state%diff_coefficient%data, current_state%local_grid, &
533  current_state%parallel%my_rank, halo_depth, involve_corners)
Here is the caller graph for this function:

◆ perform_q_diffusion()

subroutine diffusion_mod::perform_q_diffusion ( type(model_state_type), intent(inout), target  current_state,
integer, intent(in)  local_y,
integer, intent(in)  local_x 
)
private

Computes the diffusion source terms for each Q field.

Parameters
current_stateThe current model state
local_yLocal Y index
local_xLocal X index

Definition at line 433 of file diffusion.F90.

434  type(model_state_type), target, intent(inout) :: current_state
435  integer, intent(in) :: local_y, local_x
436 
437  integer :: n
438 
439  do n=1, current_state%number_q_fields
440  call general_diffusion(current_state, current_state%zq(n), current_state%sq(n), local_y, local_x, q_diffusion(:,n))
441  end do
Here is the call graph for this function:
Here is the caller graph for this function:

◆ perform_th_diffusion()

subroutine diffusion_mod::perform_th_diffusion ( type(model_state_type), intent(inout), target  current_state,
integer, intent(in)  local_y,
integer, intent(in)  local_x 
)
private

Computes the diffusion source terms for the theta field.

Parameters
current_stateThe current model state
local_yLocal Y index
local_xLocal X index

Definition at line 448 of file diffusion.F90.

449  type(model_state_type), target, intent(inout) :: current_state
450  integer, intent(in) :: local_y, local_x
451 
452  integer :: k
453  real(kind=default_precision) :: adjustment
454 
455  call general_diffusion(current_state, current_state%zth, current_state%sth, local_y, local_x, th_diffusion)
456 
457  if (current_state%use_anelastic_equations) then
458  ! This code only needs to be executed if anelastic, otherwise THREF is constant and the increment calculated here is zero
459  do k=2, current_state%local_grid%size(z_index)
460  adjustment=(current_state%global_grid%configuration%vertical%cza(k)*&
461  current_state%global_grid%configuration%vertical%dthref(k)*&
462  current_state%diff_coefficient%data(k, local_y, local_x) - current_state%global_grid%configuration%vertical%czb(k)*&
463  current_state%global_grid%configuration%vertical%dthref(k-1)*&
464  current_state%diff_coefficient%data(k-1, local_y, local_x))
465  current_state%sth%data(k, local_y, local_x)=current_state%sth%data(k, local_y, local_x)+adjustment
466  th_diffusion(k)=th_diffusion(k)+adjustment
467  end do
468  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ save_precomponent_tendencies()

subroutine diffusion_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 587 of file diffusion.F90.

588  type(model_state_type), target, intent(in) :: current_state
589  integer, intent(in) :: cxn, cyn, txn, tyn
590 
591  ! Save 3d tendency fields upon request (of 3d or profiles) and availability
592  if (l_tend_3d_th) then
593  tend_3d_th(:,tyn,txn)=current_state%sth%data(:,cyn,cxn)
594  endif
595  if (l_tend_3d_qv) then
596  tend_3d_qv(:,tyn,txn)=current_state%sq(iqv)%data(:,cyn,cxn)
597  endif
598  if (l_tend_3d_ql) then
599  tend_3d_ql(:,tyn,txn)=current_state%sq(iql)%data(:,cyn,cxn)
600  endif
601  if (l_tend_3d_qi) then
602  tend_3d_qi(:,tyn,txn)=current_state%sq(iqi)%data(:,cyn,cxn)
603  endif
604  if (l_tend_3d_qr) then
605  tend_3d_qr(:,tyn,txn)=current_state%sq(iqr)%data(:,cyn,cxn)
606  endif
607  if (l_tend_3d_qs) then
608  tend_3d_qs(:,tyn,txn)=current_state%sq(iqs)%data(:,cyn,cxn)
609  endif
610  if (l_tend_3d_qg) then
611  tend_3d_qg(:,tyn,txn)=current_state%sq(iqg)%data(:,cyn,cxn)
612  endif
613  if (l_tend_3d_tabs) then
614  tend_3d_tabs(:,tyn,txn)=current_state%sth%data(:,cyn,cxn) * current_state%global_grid%configuration%vertical%rprefrcp(:)
615  endif
616 
Here is the caller graph for this function:

◆ set_published_field_value()

subroutine diffusion_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 691 of file diffusion.F90.

692  type(component_field_value_type), intent(inout) :: field_value
693  real(kind=default_precision), dimension(:), optional :: real_1d_field
694  real(kind=default_precision), dimension(:,:), optional :: real_2d_field
695  real(kind=default_precision), dimension(:,:,:), optional :: real_3d_field
696 
697  if (present(real_1d_field)) then
698  allocate(field_value%real_1d_array(size(real_1d_field)), source=real_1d_field)
699  else if (present(real_2d_field)) then
700  allocate(field_value%real_2d_array(size(real_2d_field, 1), size(real_2d_field, 2)), source=real_2d_field)
701  else if (present(real_3d_field)) then
702  allocate(field_value%real_3d_array(size(real_3d_field, 1), size(real_3d_field, 2), size(real_3d_field, 3)), &
703  source=real_3d_field)
704  end if
Here is the caller graph for this function:

◆ timestep_callback()

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

At each timestep will compute the diffusion source terms for TH and Q fields per column if these fields are active.

Parameters
current_stateThe current model state

Definition at line 370 of file diffusion.F90.

371  type(model_state_type), target, intent(inout) :: current_state
372 
373  integer :: local_y, local_x, target_x_index, target_y_index
374 
375  local_y=current_state%column_local_y
376  local_x=current_state%column_local_x
377  target_y_index=local_y-current_state%local_grid%halo_size(y_index)
378  target_x_index=local_x-current_state%local_grid%halo_size(x_index)
379 
380  ! Zero profile tendency totals on first instance in the sum
381  if (current_state%first_timestep_column) then
382  if (l_tend_pr_tot_th) then
383  tend_pr_tot_th(:)=0.0_default_precision
384  endif
385  if (l_tend_pr_tot_qv) then
386  tend_pr_tot_qv(:)=0.0_default_precision
387  endif
388  if (l_tend_pr_tot_ql) then
389  tend_pr_tot_ql(:)=0.0_default_precision
390  endif
391  if (l_tend_pr_tot_qi) then
392  tend_pr_tot_qi(:)=0.0_default_precision
393  endif
394  if (l_tend_pr_tot_qr) then
395  tend_pr_tot_qr(:)=0.0_default_precision
396  endif
397  if (l_tend_pr_tot_qs) then
398  tend_pr_tot_qs(:)=0.0_default_precision
399  endif
400  if (l_tend_pr_tot_qg) then
401  tend_pr_tot_qg(:)=0.0_default_precision
402  endif
403  if (l_tend_pr_tot_tabs) then
404  tend_pr_tot_tabs(:)=0.0_default_precision
405  endif
406  endif ! zero totals
407 
408  if (.not. current_state%use_viscosity_and_diffusion .or. current_state%halo_column) return
409 
410  if (current_state%diffusion_halo_swap_state%swap_in_progress) then
411  ! If there is a diffusion halo swap in progress then complete it
412  call complete_nonblocking_halo_swap(current_state, current_state%diffusion_halo_swap_state, &
413  perform_local_data_copy_for_diff, copy_halo_buffer_to_diff, copy_halo_buffer_to_diff_corners)
414  end if
415 
416  if (mod(current_state%timestep, diagnostic_generation_frequency) == 0) then
417  call save_precomponent_tendencies(current_state, local_x, local_y, target_x_index, target_y_index)
418  end if
419 
420  if (current_state%th%active) call perform_th_diffusion(current_state, local_y, local_x)
421  if (current_state%number_q_fields .gt. 0) call perform_q_diffusion(current_state, local_y, local_x)
422 
423  if (mod(current_state%timestep, diagnostic_generation_frequency) == 0) then
424  call compute_component_tendencies(current_state, local_x, local_y, target_x_index, target_y_index)
425  end if
426 
Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ diagnostic_generation_frequency

integer diffusion_mod::diagnostic_generation_frequency
private

Definition at line 40 of file diffusion.F90.

40  integer :: diagnostic_generation_frequency

◆ iqg

integer diffusion_mod::iqg =0
private

Definition at line 39 of file diffusion.F90.

◆ iqi

integer diffusion_mod::iqi =0
private

Definition at line 39 of file diffusion.F90.

◆ iql

integer diffusion_mod::iql =0
private

Definition at line 39 of file diffusion.F90.

◆ iqr

integer diffusion_mod::iqr =0
private

Definition at line 39 of file diffusion.F90.

◆ iqs

integer diffusion_mod::iqs =0
private

Definition at line 39 of file diffusion.F90.

◆ iqv

integer diffusion_mod::iqv =0
private

Definition at line 39 of file diffusion.F90.

39  integer :: iqv=0, iql=0, iqr=0, iqi=0, iqs=0, iqg=0

◆ l_tend_3d_qg

logical diffusion_mod::l_tend_3d_qg
private

Definition at line 31 of file diffusion.F90.

◆ l_tend_3d_qi

logical diffusion_mod::l_tend_3d_qi
private

Definition at line 31 of file diffusion.F90.

◆ l_tend_3d_ql

logical diffusion_mod::l_tend_3d_ql
private

Definition at line 31 of file diffusion.F90.

◆ l_tend_3d_qr

logical diffusion_mod::l_tend_3d_qr
private

Definition at line 31 of file diffusion.F90.

◆ l_tend_3d_qs

logical diffusion_mod::l_tend_3d_qs
private

Definition at line 31 of file diffusion.F90.

◆ l_tend_3d_qv

logical diffusion_mod::l_tend_3d_qv
private

Definition at line 31 of file diffusion.F90.

◆ l_tend_3d_tabs

logical diffusion_mod::l_tend_3d_tabs
private

Definition at line 31 of file diffusion.F90.

◆ l_tend_3d_th

logical diffusion_mod::l_tend_3d_th
private

Definition at line 31 of file diffusion.F90.

31  logical :: l_tend_3d_th,l_tend_3d_qv,l_tend_3d_ql,l_tend_3d_qi,l_tend_3d_qr,l_tend_3d_qs,l_tend_3d_qg,l_tend_3d_tabs

◆ l_tend_pr_tot_qg

logical diffusion_mod::l_tend_pr_tot_qg
private

Definition at line 36 of file diffusion.F90.

◆ l_tend_pr_tot_qi

logical diffusion_mod::l_tend_pr_tot_qi
private

Definition at line 36 of file diffusion.F90.

◆ l_tend_pr_tot_ql

logical diffusion_mod::l_tend_pr_tot_ql
private

Definition at line 36 of file diffusion.F90.

◆ l_tend_pr_tot_qr

logical diffusion_mod::l_tend_pr_tot_qr
private

Definition at line 36 of file diffusion.F90.

◆ l_tend_pr_tot_qs

logical diffusion_mod::l_tend_pr_tot_qs
private

Definition at line 36 of file diffusion.F90.

◆ l_tend_pr_tot_qv

logical diffusion_mod::l_tend_pr_tot_qv
private

Definition at line 36 of file diffusion.F90.

◆ l_tend_pr_tot_tabs

logical diffusion_mod::l_tend_pr_tot_tabs
private

Definition at line 36 of file diffusion.F90.

◆ l_tend_pr_tot_th

logical diffusion_mod::l_tend_pr_tot_th
private

Definition at line 36 of file diffusion.F90.

36  logical :: l_tend_pr_tot_th,l_tend_pr_tot_qv,l_tend_pr_tot_ql,l_tend_pr_tot_qi,l_tend_pr_tot_qr, &
37  l_tend_pr_tot_qs,l_tend_pr_tot_qg,l_tend_pr_tot_tabs

◆ q_diffusion

real(kind=default_precision), dimension(:,:), allocatable diffusion_mod::q_diffusion
private

Definition at line 25 of file diffusion.F90.

25  real(kind=default_precision), dimension(:,:), allocatable :: q_diffusion

◆ tend_3d_qg

real(kind=default_precision), dimension(:,:,:), allocatable diffusion_mod::tend_3d_qg
private

Definition at line 29 of file diffusion.F90.

◆ tend_3d_qi

real(kind=default_precision), dimension(:,:,:), allocatable diffusion_mod::tend_3d_qi
private

Definition at line 29 of file diffusion.F90.

◆ tend_3d_ql

real(kind=default_precision), dimension(:,:,:), allocatable diffusion_mod::tend_3d_ql
private

Definition at line 29 of file diffusion.F90.

◆ tend_3d_qr

real(kind=default_precision), dimension(:,:,:), allocatable diffusion_mod::tend_3d_qr
private

Definition at line 29 of file diffusion.F90.

◆ tend_3d_qs

real(kind=default_precision), dimension(:,:,:), allocatable diffusion_mod::tend_3d_qs
private

Definition at line 29 of file diffusion.F90.

◆ tend_3d_qv

real(kind=default_precision), dimension(:,:,:), allocatable diffusion_mod::tend_3d_qv
private

Definition at line 29 of file diffusion.F90.

◆ tend_3d_tabs

real(kind=default_precision), dimension(:,:,:), allocatable diffusion_mod::tend_3d_tabs
private

Definition at line 29 of file diffusion.F90.

◆ tend_3d_th

real(kind=default_precision), dimension(:,:,:), allocatable diffusion_mod::tend_3d_th
private

Definition at line 29 of file diffusion.F90.

29  real(kind=default_precision), dimension(:,:,:), allocatable :: &
30  tend_3d_th,tend_3d_qv,tend_3d_ql,tend_3d_qi,tend_3d_qr,tend_3d_qs,tend_3d_qg,tend_3d_tabs

◆ tend_pr_tot_qg

real(kind=default_precision), dimension(:), allocatable diffusion_mod::tend_pr_tot_qg
private

Definition at line 33 of file diffusion.F90.

◆ tend_pr_tot_qi

real(kind=default_precision), dimension(:), allocatable diffusion_mod::tend_pr_tot_qi
private

Definition at line 33 of file diffusion.F90.

◆ tend_pr_tot_ql

real(kind=default_precision), dimension(:), allocatable diffusion_mod::tend_pr_tot_ql
private

Definition at line 33 of file diffusion.F90.

◆ tend_pr_tot_qr

real(kind=default_precision), dimension(:), allocatable diffusion_mod::tend_pr_tot_qr
private

Definition at line 33 of file diffusion.F90.

◆ tend_pr_tot_qs

real(kind=default_precision), dimension(:), allocatable diffusion_mod::tend_pr_tot_qs
private

Definition at line 33 of file diffusion.F90.

◆ tend_pr_tot_qv

real(kind=default_precision), dimension(:), allocatable diffusion_mod::tend_pr_tot_qv
private

Definition at line 33 of file diffusion.F90.

◆ tend_pr_tot_tabs

real(kind=default_precision), dimension(:), allocatable diffusion_mod::tend_pr_tot_tabs
private

Definition at line 33 of file diffusion.F90.

◆ tend_pr_tot_th

real(kind=default_precision), dimension(:), allocatable diffusion_mod::tend_pr_tot_th
private

Definition at line 33 of file diffusion.F90.

33  real(kind=default_precision), dimension(:), allocatable :: &
34  tend_pr_tot_th,tend_pr_tot_qv,tend_pr_tot_ql,tend_pr_tot_qi,tend_pr_tot_qr,tend_pr_tot_qs,tend_pr_tot_qg, &
35  tend_pr_tot_tabs

◆ th_diffusion

real(kind=default_precision), dimension(:), allocatable diffusion_mod::th_diffusion
private

Definition at line 24 of file diffusion.F90.

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