MONC
def_tvd_diagnostic_terms.F90
Go to the documentation of this file.
2 
3  use state_mod, only : model_state_type
6 
7  implicit none
8 
10 
11  ! TVD advective fluxes of u,v, w, q through the bottom face of grid
12  ! Needed in profile_diagnostics to calc the turbulent flux diagnostic
13  ! when TVD advection is used
14  real(kind=default_precision), dimension(:,:,:), allocatable :: adv_u_dgs
15  real(kind=default_precision), dimension(:,:,:), allocatable :: adv_v_dgs
16  real(kind=default_precision), dimension(:,:,:), allocatable :: adv_w_dgs
17  real(kind=default_precision), dimension(:,:,:), allocatable :: adv_th_dgs
18 
19  real(kind=default_precision), dimension(:,:,:,:), allocatable :: adv_q_dgs
20 
22 
24 
25 contains
26 
27  subroutine allocate_tvd_diagnostic_terms(current_state, tvd_dgs_terms)
28 
29  implicit none
30 
31  type(model_state_type), target, intent(in) :: current_state
32  type(str_tvd_diagnostic_terms), intent(inout) :: tvd_dgs_terms
33 
34  integer :: k_top, x_local, y_local
35 
36  k_top = current_state%local_grid%size(z_index) + current_state%local_grid%halo_size(z_index) * 2
37  x_local = current_state%local_grid%size(x_index) + current_state%local_grid%halo_size(x_index) * 2
38  y_local = current_state%local_grid%size(y_index) + current_state%local_grid%halo_size(x_index) * 2
39 
40  allocate(tvd_dgs_terms%adv_u_dgs(k_top, y_local, x_local), &
41  tvd_dgs_terms%adv_v_dgs(k_top, y_local, x_local), &
42  tvd_dgs_terms%adv_w_dgs(k_top, y_local, x_local), &
43  tvd_dgs_terms%adv_th_dgs(k_top, y_local, x_local))
44 
45  tvd_dgs_terms%adv_u_dgs(:,:,:) = 0.0
46  tvd_dgs_terms%adv_v_dgs(:,:,:) = 0.0
47  tvd_dgs_terms%adv_w_dgs(:,:,:) = 0.0
48  tvd_dgs_terms%adv_th_dgs(:,:,:) = 0.0
49 
50  if (current_state%number_q_fields > 0) then
51  allocate(tvd_dgs_terms%adv_q_dgs(k_top, y_local, x_local, current_state%number_q_fields))
52  tvd_dgs_terms%adv_q_dgs(:,:,:,:)= 0.0
53  endif
54 
55  end subroutine allocate_tvd_diagnostic_terms
56 
57  subroutine deallocate_tvd_diagnostic_terms(current_state, tvd_dgs_terms)
58  type(model_state_type), target, intent(in) :: current_state
59  type(str_tvd_diagnostic_terms), intent(inout) :: tvd_dgs_terms
60 
61  deallocate(tvd_dgs_terms%adv_u_dgs, &
62  tvd_dgs_terms%adv_v_dgs, &
63  tvd_dgs_terms%adv_w_dgs, &
64  tvd_dgs_terms%adv_th_dgs)
65 
66  if (current_state%number_q_fields > 0) then
67  deallocate(tvd_dgs_terms%adv_q_dgs)
68  endif
69 
70  end subroutine deallocate_tvd_diagnostic_terms
71 
72 
73 end module def_tvd_diagnostic_terms
def_tvd_diagnostic_terms::str_tvd_diagnostic_terms
Definition: def_tvd_diagnostic_terms.F90:9
grids_mod::x_index
integer, parameter, public x_index
Definition: grids.F90:14
grids_mod::y_index
integer, parameter, public y_index
Definition: grids.F90:14
grids_mod::global_grid_type
Defines the global grid.
Definition: grids.F90:107
def_tvd_diagnostic_terms
Definition: def_tvd_diagnostic_terms.F90:1
grids_mod::z_index
integer, parameter, public z_index
Grid index parameters.
Definition: grids.F90:14
def_tvd_diagnostic_terms::tvd_dgs_terms
type(str_tvd_diagnostic_terms) tvd_dgs_terms
Definition: def_tvd_diagnostic_terms.F90:23
state_mod::model_state_type
The ModelState which represents the current state of a run.
Definition: state.F90:39
grids_mod::local_grid_type
Defined the local grid, i.e. the grid held on this process after decomposition.
Definition: grids.F90:118
datadefn_mod
Contains common definitions for the data and datatypes used by MONC.
Definition: datadefn.F90:2
def_tvd_diagnostic_terms::allocate_tvd_diagnostic_terms
subroutine allocate_tvd_diagnostic_terms(current_state, tvd_dgs_terms)
Definition: def_tvd_diagnostic_terms.F90:28
grids_mod
Functionality to support the different types of grid and abstraction between global grids and local o...
Definition: grids.F90:5
def_tvd_diagnostic_terms::deallocate_tvd_diagnostic_terms
subroutine deallocate_tvd_diagnostic_terms(current_state, tvd_dgs_terms)
Definition: def_tvd_diagnostic_terms.F90:58
datadefn_mod::default_precision
integer, parameter, public default_precision
MPI communication type which we use for the prognostic and calculation data.
Definition: datadefn.F90:17
state_mod
The model state which represents the current state of a run.
Definition: state.F90:2