MONC
Functions/Subroutines | Variables
thadvection_mod Module Reference

Specific theta advection, which involves the vertical advection of reference state and advection of mean baroclinicity. More...

Functions/Subroutines

type(component_descriptor_type) function, public thadvection_get_descriptor ()
 Provides a description of this component for the core to register. More...
 
subroutine initialisation_callback (current_state)
 Initialisation callback to set up the variables and data needed by the component. More...
 
subroutine finalisation_callback (current_state)
 
subroutine timestep_callback (current_state)
 Timestep callback, will call the two separate procedures to do their advection if needed. More...
 
subroutine vertical_advection_of_reference_state (current_state, local_y, local_x)
 Vertical advection of the reference state. It doesn't seem consistent to do the advection in this way if TVD advection of the deviation from the reference state has been selected. Separate vertical advection of the reference state was introduced to improve energy conservation when carrying out idealized gravity wave simulations in a deep, dry isothermal layer, for which the difference in potential temp between top and bottom was of order 100K. In less extreme cases the benefits are unlikely to be significant and with TVD advection energy conservation has been compromised so the best way forward might be to recombine the reference state into l_th. More...
 
subroutine advection_of_mean_baroclinicity (current_state, local_y, local_x)
 Performs advection of the mean baroclinicity if appropriate. 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 field_information_retrieval_callback (current_state, name, field_information)
 Field information retrieval callback, this returns information for a specific component's 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 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

logical baroclinicity_use_geostrophic_shear
 
logical l_advect_mean_baroclinicity
 
real(kind=default_precision) fcoriol
 
real(kind=default_precision) fcoriol_over_g
 
real(kind=default_precision) rate_change_geostrophic_wind_x
 
real(kind=default_precision) rate_change_geostrophic_wind_y
 
real(kind=default_precision) multiplicative_factor_x
 
real(kind=default_precision) multiplicative_factor_y
 
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_th
 
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_tabs
 
logical l_tend_3d_th
 
logical l_tend_3d_tabs
 
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_th
 
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_tabs
 
logical l_tend_pr_tot_th
 
logical l_tend_pr_tot_tabs
 
integer diagnostic_generation_frequency
 

Detailed Description

Specific theta advection, which involves the vertical advection of reference state and advection of mean baroclinicity.

Function/Subroutine Documentation

◆ advection_of_mean_baroclinicity()

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

Performs advection of the mean baroclinicity if appropriate.

Parameters
current_stateThe current model state
local_yThe local Y of the column
local_xThe local X of the column

Definition at line 201 of file thadvection.F90.

202  type(model_state_type), target, intent(inout) :: current_state
203  integer, intent(in) :: local_y, local_x
204 
205  integer :: k
206 
207  if (l_advect_mean_baroclinicity) then
208  if (current_state%use_anelastic_equations) then
209  do k=2, current_state%local_grid%size(z_index)
210  current_state%sth%data(k, local_y, local_x)=current_state%sth%data(k, local_y, local_x)+&
211  current_state%global_grid%configuration%vertical%thref(k)*fcoriol_over_g*&
212  ((current_state%v%data(k, local_y, local_x) + current_state%vgal) * rate_change_geostrophic_wind_x-&
213  (current_state%u%data(k, local_y, local_x) + current_state%ugal) * rate_change_geostrophic_wind_y)
214  end do
215  else
216  do k=2, current_state%local_grid%size(z_index)
217  current_state%sth%data(k, local_y, local_x)=current_state%sth%data(k, local_y, local_x)+&
218  ((current_state%v%data(k, local_y, local_x) + current_state%vgal) * multiplicative_factor_x-&
219  (current_state%u%data(k, local_y, local_x) + current_state%ugal) * multiplicative_factor_y)
220  end do
221  end if
222  end if
223 
Here is the caller graph for this function:

◆ compute_component_tendencies()

subroutine thadvection_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 254 of file thadvection.F90.

255  type(model_state_type), target, intent(inout) :: current_state
256  integer, intent(in) :: cxn, cyn, txn, tyn
257 
258  ! Calculate change in tendency due to component
259  if (l_tend_3d_th) then
260  tend_3d_th(:,tyn,txn)=current_state%sth%data(:,cyn,cxn) - tend_3d_th(:,tyn,txn)
261  endif
262  if (l_tend_3d_tabs) then
263  tend_3d_tabs(:,tyn,txn)= &
264  current_state%sth%data(:,cyn,cxn) * current_state%global_grid%configuration%vertical%rprefrcp(:) &
265  - tend_3d_tabs(:,tyn,txn)
266  endif
267 
268  ! Add local tendency fields to the profile total
269  if (l_tend_pr_tot_th) then
270  tend_pr_tot_th(:)=tend_pr_tot_th(:) + tend_3d_th(:,tyn,txn)
271  endif
272  if (l_tend_pr_tot_tabs) then
273  tend_pr_tot_tabs(:)=tend_pr_tot_tabs(:) + tend_3d_tabs(:,tyn,txn)
274  endif
275 
Here is the caller graph for this function:

◆ field_information_retrieval_callback()

subroutine thadvection_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 component's published field.

Parameters
current_stateCurrent model state
nameThe name of the field to retrieve information for
field_informationPopulated with information about the field
strcompStarting index within 1st argument string that matches substring (2nd argument); 0 if not a match.

Definition at line 284 of file thadvection.F90.

285  type(model_state_type), target, intent(inout) :: current_state
286  character(len=*), intent(in) :: name
287  type(component_field_information_type), intent(out) :: field_information
288  integer :: strcomp
289 
290  ! Field information for 3d
291  strcomp=index(name, "thadvection_3d_local")
292  if (strcomp .ne. 0) then
293  field_information%field_type=component_array_field_type
294  field_information%number_dimensions=3
295  field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
296  field_information%dimension_sizes(2)=current_state%local_grid%size(y_index)
297  field_information%dimension_sizes(3)=current_state%local_grid%size(x_index)
298  field_information%data_type=component_double_data_type
299 
300  if (name .eq. "tend_th_thadvection_3d_local") then
301  field_information%enabled=l_tend_3d_th
302  else if (name .eq. "tend_tabs_thadvection_3d_local") then
303  field_information%enabled=l_tend_3d_tabs
304  else
305  field_information%enabled=.true.
306  end if
307 
308  end if !end 3d check
309 
310  ! Field information for profiles
311  strcomp=index(name, "thadvection_profile_total_local")
312  if (strcomp .ne. 0) then
313  field_information%field_type=component_array_field_type
314  field_information%number_dimensions=1
315  field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
316  field_information%data_type=component_double_data_type
317 
318  if (name .eq. "tend_th_thadvection_profile_total_local") then
319  field_information%enabled=l_tend_pr_tot_th
320  else if (name .eq. "tend_tabs_thadvection_profile_total_local") then
321  field_information%enabled=l_tend_pr_tot_tabs
322  else
323  field_information%enabled=.true.
324  end if
325 
326  end if !end profile check
327 
Here is the caller graph for this function:

◆ field_value_retrieval_callback()

subroutine thadvection_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 335 of file thadvection.F90.

336  type(model_state_type), target, intent(inout) :: current_state
337  character(len=*), intent(in) :: name
338  type(component_field_value_type), intent(out) :: field_value
339 
340  ! 3d Tendency Fields
341  if (name .eq. "tend_th_thadvection_3d_local" .and. allocated(tend_3d_th)) then
342  call set_published_field_value(field_value, real_3d_field=tend_3d_th)
343  else if (name .eq. "tend_tabs_thadvection_3d_local" .and. allocated(tend_3d_tabs)) then
344  call set_published_field_value(field_value, real_3d_field=tend_3d_tabs)
345 
346  ! Profile Tendency Fields
347  else if (name .eq. "tend_th_thadvection_profile_total_local" .and. allocated(tend_pr_tot_th)) then
348  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_th)
349  else if (name .eq. "tend_tabs_thadvection_profile_total_local" .and. allocated(tend_pr_tot_tabs)) then
350  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_tabs)
351  end if
352 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ finalisation_callback()

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

Definition at line 119 of file thadvection.F90.

120  type(model_state_type), target, intent(inout) :: current_state
121 
122  if (allocated(tend_3d_th)) deallocate(tend_3d_th)
123  if (allocated(tend_3d_tabs)) deallocate(tend_3d_tabs)
124 
125  if (allocated(tend_pr_tot_th)) deallocate(tend_pr_tot_th)
126  if (allocated(tend_pr_tot_tabs)) deallocate(tend_pr_tot_tabs)
127 
Here is the caller graph for this function:

◆ initialisation_callback()

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

Initialisation callback to set up the variables and data needed by the component.

Parameters
current_stateThe current model state

Definition at line 60 of file thadvection.F90.

61  type(model_state_type), target, intent(inout) :: current_state
62 
63  !AH (08/02/18) - by default l_advect_mean_baroclinicity is false. It
64  ! is switch to true if passive_q true and baroclinicity_use_geostrophic_shear is true
65  ! If passive q is false code will run but a warning is issued.
66  ! Mean advection of baroclinicity only works for dry case
67  baroclinicity_use_geostrophic_shear=options_get_logical(current_state%options_database, "baroclinicity_use_geostrophic_shear")
68  if (baroclinicity_use_geostrophic_shear) then
69  if (current_state%passive_q) then
70  l_advect_mean_baroclinicity = .true.
71  else
72  call log_master_log(log_warn, "The combination if baroclinicity and active q is not allowed, code will run but"// &
73  " no advection of mean baroclinicity")
74  l_advect_mean_baroclinicity = .false.
75  endif
76  endif
77  fcoriol=options_get_real(current_state%options_database, "fcoriol")
78  rate_change_geostrophic_wind_x=options_get_real(current_state%options_database, "rate_change_geostrophic_wind_x")
79  rate_change_geostrophic_wind_y=options_get_real(current_state%options_database, "rate_change_geostrophic_wind_y")
80  fcoriol_over_g = fcoriol/g
81  multiplicative_factor_x=rate_change_geostrophic_wind_x*current_state%thref0*fcoriol_over_g
82  multiplicative_factor_y=rate_change_geostrophic_wind_y*current_state%thref0*fcoriol_over_g
83 
84  ! Set tendency diagnostic logicals based on availability
85  ! Need to use 3d tendencies to compute the profiles, so they will be allocated
86  ! in the case where profiles are available
87  l_tend_pr_tot_th = current_state%th%active
88  l_tend_pr_tot_tabs = l_tend_pr_tot_th
89 
90  l_tend_3d_th = current_state%th%active .or. l_tend_pr_tot_th
91  l_tend_3d_tabs = l_tend_3d_th
92 
93  ! Allocate 3d tendency fields upon availability
94  if (l_tend_3d_th) then
95  allocate( tend_3d_th(current_state%local_grid%size(z_index), &
96  current_state%local_grid%size(y_index), &
97  current_state%local_grid%size(x_index) ) )
98  endif
99  if (l_tend_3d_tabs) then
100  allocate( tend_3d_tabs(current_state%local_grid%size(z_index), &
101  current_state%local_grid%size(y_index), &
102  current_state%local_grid%size(x_index) ) )
103  endif
104 
105  ! Allocate profile tendency fields upon availability
106  if (l_tend_pr_tot_th) then
107  allocate( tend_pr_tot_th(current_state%local_grid%size(z_index)) )
108  endif
109  if (l_tend_pr_tot_tabs) then
110  allocate( tend_pr_tot_tabs(current_state%local_grid%size(z_index)) )
111  endif
112 
113  ! Save the sampling_frequency to force diagnostic calculation on select time steps
114  diagnostic_generation_frequency=options_get_integer(current_state%options_database, "sampling_frequency")
115 
Here is the caller graph for this function:

◆ save_precomponent_tendencies()

subroutine thadvection_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 233 of file thadvection.F90.

234  type(model_state_type), target, intent(in) :: current_state
235  integer, intent(in) :: cxn, cyn, txn, tyn
236 
237  ! Save 3d tendency fields upon request (of 3d or profiles) and availability
238  if (l_tend_3d_th) then
239  tend_3d_th(:,tyn,txn)=current_state%sth%data(:,cyn,cxn)
240  endif
241  if (l_tend_3d_tabs) then
242  tend_3d_tabs(:,tyn,txn)=current_state%sth%data(:,cyn,cxn) * current_state%global_grid%configuration%vertical%rprefrcp(:)
243  endif
244 
Here is the caller graph for this function:

◆ set_published_field_value()

subroutine thadvection_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 360 of file thadvection.F90.

361  type(component_field_value_type), intent(inout) :: field_value
362  real(kind=default_precision), dimension(:), optional :: real_1d_field
363  real(kind=default_precision), dimension(:,:), optional :: real_2d_field
364  real(kind=default_precision), dimension(:,:,:), optional :: real_3d_field
365 
366  if (present(real_1d_field)) then
367  allocate(field_value%real_1d_array(size(real_1d_field)), source=real_1d_field)
368  else if (present(real_2d_field)) then
369  allocate(field_value%real_2d_array(size(real_2d_field, 1), size(real_2d_field, 2)), source=real_2d_field)
370  else if (present(real_3d_field)) then
371  allocate(field_value%real_3d_array(size(real_3d_field, 1), size(real_3d_field, 2), size(real_3d_field, 3)), &
372  source=real_3d_field)
373  end if
Here is the caller graph for this function:

◆ thadvection_get_descriptor()

type(component_descriptor_type) function, public thadvection_mod::thadvection_get_descriptor

Provides a description of this component for the core to register.

Returns
The descriptor containing registration information for this component

Definition at line 39 of file thadvection.F90.

40  thadvection_get_descriptor%name="th_advection"
41  thadvection_get_descriptor%version=0.1
42  thadvection_get_descriptor%initialisation=>initialisation_callback
43  thadvection_get_descriptor%timestep=>timestep_callback
44  thadvection_get_descriptor%finalisation=>finalisation_callback
45 
46  thadvection_get_descriptor%field_value_retrieval=>field_value_retrieval_callback
47  thadvection_get_descriptor%field_information_retrieval=>field_information_retrieval_callback
48  allocate(thadvection_get_descriptor%published_fields(2+2))
49 
50  thadvection_get_descriptor%published_fields(1)= "tend_th_thadvection_3d_local"
51  thadvection_get_descriptor%published_fields(2)= "tend_tabs_thadvection_3d_local"
52 
53  thadvection_get_descriptor%published_fields(2+1)= "tend_th_thadvection_profile_total_local"
54  thadvection_get_descriptor%published_fields(2+2)= "tend_tabs_thadvection_profile_total_local"
55 
Here is the call graph for this function:

◆ timestep_callback()

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

Timestep callback, will call the two separate procedures to do their advection if needed.

Parameters
current_stateThe current model state
target_(x/y)_indexThis is the index with the halos subtracted. This is needed so that diagnostic does not include halos and to prevent array out-of-bounds

Definition at line 135 of file thadvection.F90.

136  type(model_state_type), target, intent(inout) :: current_state
137 
138  integer :: current_x_index, current_y_index, target_x_index, target_y_index
139 
140  current_x_index=current_state%column_local_x
141  current_y_index=current_state%column_local_y
142  target_y_index=current_y_index-current_state%local_grid%halo_size(y_index)
143  target_x_index=current_x_index-current_state%local_grid%halo_size(x_index)
144 
145  ! Zero profile tendency totals on first instance in the sum
146  if (current_state%first_timestep_column) then
147  if (l_tend_pr_tot_th) then
148  tend_pr_tot_th(:)=0.0_default_precision
149  endif
150  if (l_tend_pr_tot_tabs) then
151  tend_pr_tot_tabs(:)=0.0_default_precision
152  endif
153  endif ! zero totals
154 
155  if (mod(current_state%timestep, diagnostic_generation_frequency) == 0 .and. .not. current_state%halo_column) then
156  call save_precomponent_tendencies(current_state, current_x_index, current_y_index, target_x_index, target_y_index)
157  end if
158 
159  call vertical_advection_of_reference_state(current_state, current_state%column_local_y, current_state%column_local_x)
160  call advection_of_mean_baroclinicity(current_state, current_state%column_local_y, current_state%column_local_x)
161 
162  if (mod(current_state%timestep, diagnostic_generation_frequency) == 0 .and. .not. current_state%halo_column) then
163  call compute_component_tendencies(current_state, current_x_index, current_y_index, target_x_index, target_y_index)
164  end if
165 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ vertical_advection_of_reference_state()

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

Vertical advection of the reference state. It doesn't seem consistent to do the advection in this way if TVD advection of the deviation from the reference state has been selected. Separate vertical advection of the reference state was introduced to improve energy conservation when carrying out idealized gravity wave simulations in a deep, dry isothermal layer, for which the difference in potential temp between top and bottom was of order 100K. In less extreme cases the benefits are unlikely to be significant and with TVD advection energy conservation has been compromised so the best way forward might be to recombine the reference state into l_th.

Parameters
current_stateThe current model state
local_yThe local Y of the column
local_xThe local X of the column

Definition at line 177 of file thadvection.F90.

178  type(model_state_type), target, intent(inout) :: current_state
179  integer, intent(in) :: local_y, local_x
180 
181  integer :: k
182  real(kind=default_precision) :: sctmp1, sctmp2
183 
184  if (current_state%use_anelastic_equations) then
185  ! This code only needs to be executed if anelastic, otherwise THREF is constant and the increment calculated here is zero
186  do k=2, current_state%local_grid%size(z_index)
187  sctmp1=current_state%global_grid%configuration%vertical%tzc1(k)*2.0_default_precision*&
188  current_state%global_grid%configuration%vertical%dthref(k-1)
189  sctmp2=current_state%global_grid%configuration%vertical%tzc2(k)*2.0_default_precision*&
190  current_state%global_grid%configuration%vertical%dthref(k)
191  current_state%sth%data(k, local_y, local_x)=current_state%sth%data(k, local_y, local_x)-(sctmp1*&
192  current_state%w%data(k-1, local_y, local_x) + sctmp2*current_state%w%data(k, local_y, local_x))
193  end do
194  end if
Here is the caller graph for this function:

Variable Documentation

◆ baroclinicity_use_geostrophic_shear

logical thadvection_mod::baroclinicity_use_geostrophic_shear
private

Definition at line 18 of file thadvection.F90.

18  logical :: baroclinicity_use_geostrophic_shear

◆ diagnostic_generation_frequency

integer thadvection_mod::diagnostic_generation_frequency
private

Definition at line 31 of file thadvection.F90.

31  integer :: diagnostic_generation_frequency

◆ fcoriol

real(kind=default_precision) thadvection_mod::fcoriol
private

Definition at line 20 of file thadvection.F90.

20  real(kind=default_precision) :: fcoriol, fcoriol_over_g, rate_change_geostrophic_wind_x, rate_change_geostrophic_wind_y, &
21  multiplicative_factor_x, multiplicative_factor_y

◆ fcoriol_over_g

real(kind=default_precision) thadvection_mod::fcoriol_over_g
private

Definition at line 20 of file thadvection.F90.

◆ l_advect_mean_baroclinicity

logical thadvection_mod::l_advect_mean_baroclinicity
private

Definition at line 19 of file thadvection.F90.

19  logical :: l_advect_mean_baroclinicity

◆ l_tend_3d_tabs

logical thadvection_mod::l_tend_3d_tabs
private

Definition at line 26 of file thadvection.F90.

◆ l_tend_3d_th

logical thadvection_mod::l_tend_3d_th
private

Definition at line 26 of file thadvection.F90.

26  logical :: l_tend_3d_th, l_tend_3d_tabs

◆ l_tend_pr_tot_tabs

logical thadvection_mod::l_tend_pr_tot_tabs
private

Definition at line 29 of file thadvection.F90.

◆ l_tend_pr_tot_th

logical thadvection_mod::l_tend_pr_tot_th
private

Definition at line 29 of file thadvection.F90.

29  logical :: l_tend_pr_tot_th, l_tend_pr_tot_tabs

◆ multiplicative_factor_x

real(kind=default_precision) thadvection_mod::multiplicative_factor_x
private

Definition at line 20 of file thadvection.F90.

◆ multiplicative_factor_y

real(kind=default_precision) thadvection_mod::multiplicative_factor_y
private

Definition at line 20 of file thadvection.F90.

◆ rate_change_geostrophic_wind_x

real(kind=default_precision) thadvection_mod::rate_change_geostrophic_wind_x
private

Definition at line 20 of file thadvection.F90.

◆ rate_change_geostrophic_wind_y

real(kind=default_precision) thadvection_mod::rate_change_geostrophic_wind_y
private

Definition at line 20 of file thadvection.F90.

◆ tend_3d_tabs

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

Definition at line 25 of file thadvection.F90.

◆ tend_3d_th

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

Definition at line 25 of file thadvection.F90.

25  real(kind=default_precision), dimension(:,:,:), allocatable :: tend_3d_th, tend_3d_tabs

◆ tend_pr_tot_tabs

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

Definition at line 28 of file thadvection.F90.

◆ tend_pr_tot_th

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

Definition at line 28 of file thadvection.F90.

28  real(kind=default_precision), dimension(:), allocatable :: tend_pr_tot_th, tend_pr_tot_tabs
logging_mod::log_warn
integer, parameter, public log_warn
Log WARNING and ERROR messages.
Definition: logging.F90:12
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
optionsdatabase_mod::options_get_logical
logical function, public options_get_logical(options_database, key, index)
Retrieves a logical value from the database that matches the provided key.
Definition: optionsdatabase.F90:154
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