MONC
Functions/Subroutines | Variables
diagnostics_mod Module Reference

Functions/Subroutines

type(component_descriptor_type) function, public diagnostics_get_descriptor ()
 Provides the component descriptor for the core to register. More...
 
subroutine initialisation_callback (current_state)
 
subroutine timestep_callback (current_state)
 
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...
 

Variables

integer total_points
 
integer hqlmax
 
integer ncl_col
 
real(kind=default_precision), dimension(:), allocatable tempfac
 
real(kind=default_precision), dimension(:), allocatable theta_tot
 
real(kind=default_precision) totqv
 
real(kind=default_precision) totql
 
real(kind=default_precision) totqz
 
real(kind=default_precision) wmax
 
real(kind=default_precision) wmin
 
real(kind=default_precision) qlmax
 
real(kind=default_precision) qlcrit
 
real(kind=default_precision) cltop_av
 
real(kind=default_precision) clbas_av
 
real(kind=default_precision) cltop
 
real(kind=default_precision) clbas
 

Function/Subroutine Documentation

◆ diagnostics_get_descriptor()

type(component_descriptor_type) function, public diagnostics_mod::diagnostics_get_descriptor ( )

Provides the component descriptor for the core to register.

Returns
The descriptor describing this component

Definition at line 27 of file diagnostics.F90.

27  diagnostics_get_descriptor%name="diagnostics"
28  diagnostics_get_descriptor%version=0.1
29 
30  diagnostics_get_descriptor%initialisation=>initialisation_callback
31  diagnostics_get_descriptor%timestep=>timestep_callback
32 
33  diagnostics_get_descriptor%field_value_retrieval=>field_value_retrieval_callback
34  diagnostics_get_descriptor%field_information_retrieval=>field_information_retrieval_callback
35  allocate(diagnostics_get_descriptor%published_fields(15))
36 
37  diagnostics_get_descriptor%published_fields(1)="totpnts"
38  diagnostics_get_descriptor%published_fields(2)="zn_top"
39  diagnostics_get_descriptor%published_fields(3)="wmax_local"
40  diagnostics_get_descriptor%published_fields(4)="wmin_local"
41  diagnostics_get_descriptor%published_fields(5)="totqv_local"
42  diagnostics_get_descriptor%published_fields(6)="totql_local"
43  diagnostics_get_descriptor%published_fields(7)="totqz_local"
44  diagnostics_get_descriptor%published_fields(8)="qlmax_local"
45  diagnostics_get_descriptor%published_fields(9)="hqlmax_local"
46  diagnostics_get_descriptor%published_fields(10)="cltop_local"
47  diagnostics_get_descriptor%published_fields(11)="clbas_local"
48  diagnostics_get_descriptor%published_fields(12)="cltop_av_local"
49  diagnostics_get_descriptor%published_fields(13)="clbas_av_local"
50  diagnostics_get_descriptor%published_fields(14)="ncl_col_local"
51  diagnostics_get_descriptor%published_fields(15)="thdiag_local"
Here is the call graph for this function:

◆ field_information_retrieval_callback()

subroutine diagnostics_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 145 of file diagnostics.F90.

145  type(model_state_type), target, intent(inout) :: current_state
146  character(len=*), intent(in) :: name
147  type(component_field_information_type), intent(out) :: field_information
148 
149  field_information%field_type=component_scalar_field_type
150  field_information%data_type=component_double_data_type
151  field_information%enabled=.true.
152  if (name .eq. "thdiag_local") then
153  field_information%field_type=component_array_field_type
154  field_information%number_dimensions=1
155  field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
156  endif
157 
158  if (name .eq. "totpnts" .or. name .eq. "hqlmax_local" .or. name .eq. "ncl_col_local") then
159  field_information%data_type=component_integer_data_type
160  else
161  field_information%data_type=component_double_data_type
162  end if
Here is the caller graph for this function:

◆ field_value_retrieval_callback()

subroutine diagnostics_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 170 of file diagnostics.F90.

170  type(model_state_type), target, intent(inout) :: current_state
171  character(len=*), intent(in) :: name
172  type(component_field_value_type), intent(out) :: field_value
173 
174  integer :: i
175 
176  if (name .eq. "totpnts") then
177  field_value%scalar_int=total_points
178  else if (name .eq. "zn_top") then
179  field_value%scalar_real=current_state%global_grid%configuration%vertical%zn(current_state%global_grid%size(z_index))
180  else if (name .eq. "wmax_local") then
181  field_value%scalar_real=wmax
182  else if (name .eq. "wmin_local") then
183  field_value%scalar_real=wmin
184  else if (name .eq. "totqv_local") then
185  field_value%scalar_real=totqv
186  else if (name .eq. "totql_local") then
187  field_value%scalar_real=totql
188  else if (name .eq. "totqz_local") then
189  field_value%scalar_real=totqz
190  else if (name .eq. "qlmax_local") then
191  field_value%scalar_real=qlmax
192  else if (name .eq. "hqlmax_local") then
193  field_value%scalar_int=hqlmax
194  else if (name .eq. "cltop_local") then
195  field_value%scalar_real=cltop
196  else if (name .eq. "clbas_local") then
197  field_value%scalar_real=clbas
198  else if (name .eq. "cltop_av_local") then
199  field_value%scalar_real=cltop_av
200  else if (name .eq. "clbas_av_local") then
201  field_value%scalar_real=clbas_av
202  else if (name .eq. "ncl_col_local") then
203  field_value%scalar_int=ncl_col
204  else if (name .eq. "thdiag_local") then
205  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
206  do i=1, current_state%local_grid%size(z_index)
207  field_value%real_1d_array(i)=maxval(current_state%th%data(i,:,:) + current_state%thref0)
208  enddo
209  end if
Here is the caller graph for this function:

◆ initialisation_callback()

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

Definition at line 55 of file diagnostics.F90.

55  type(model_state_type), target, intent(inout) :: current_state
56 
57  integer :: k
58 
59  qlcrit=options_get_real(current_state%options_database, "qlcrit")
60 
61  total_points=current_state%global_grid%size(z_index) * current_state%global_grid%size(y_index) * &
62  current_state%global_grid%size(x_index)
63 
64  allocate(tempfac(current_state%local_grid%size(z_index)), theta_tot(current_state%local_grid%size(z_index)))
65  do k=2, current_state%local_grid%size(z_index)
66  tempfac(k)=current_state%global_grid%configuration%vertical%dz(k)*&
67  current_state%global_grid%configuration%vertical%rhon(k)/total_points
68  end do
Here is the caller graph for this function:

◆ timestep_callback()

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

Definition at line 72 of file diagnostics.F90.

72  type(model_state_type), target, intent(inout) :: current_state
73 
74  integer :: k, i
75  real(kind=DEFAULT_PRECISION) :: cltop_col, clbas_col
76 
77  if (current_state%first_timestep_column) then
78  totqv=0.0_default_precision
79  totql=0.0_default_precision
80  totqz=0.0_default_precision
81  wmax=0.0_default_precision
82  wmin=0.0_default_precision
83  qlmax=0.0_default_precision
84  cltop_av=0.0_default_precision
85  clbas_av=0.0_default_precision
86  cltop=0.0_default_precision
87  clbas=0.0_default_precision
88  ncl_col=0
89  end if
90  if (.not. current_state%halo_column) then
91  cltop_col=0.0_default_precision
92  clbas_col=0.0_default_precision
93  do k=2, current_state%local_grid%size(z_index)
94  wmax=max(wmax, current_state%w%data(k, current_state%column_local_y, current_state%column_local_x))
95  wmin=min(wmin, current_state%w%data(k, current_state%column_local_y, current_state%column_local_x))
96 
97  if (current_state%number_q_fields .gt. 0) then
98  if (current_state%liquid_water_mixing_ratio_index .gt. 0 .and. &
99  current_state%number_q_fields .ge. current_state%liquid_water_mixing_ratio_index) then
100  if (qlmax .lt. current_state%q(current_state%liquid_water_mixing_ratio_index)%data(k, &
101  current_state%column_local_y, current_state%column_local_x)) then
102  qlmax=max(qlmax, current_state%q(current_state%liquid_water_mixing_ratio_index)%data(k, &
103  current_state%column_local_y, current_state%column_local_x))
104  hqlmax=k
105  end if
106 
107  if (current_state%q(current_state%liquid_water_mixing_ratio_index)%data(k, &
108  current_state%column_local_y, current_state%column_local_x) .gt. qlcrit) then
109  cltop_col=current_state%global_grid%configuration%vertical%zn(k)
110  cltop=max(cltop, current_state%global_grid%configuration%vertical%zn(k))
111  clbas=min(clbas, current_state%global_grid%configuration%vertical%zn(k))
112  end if
113 
114  if (current_state%q(current_state%liquid_water_mixing_ratio_index)%data(current_state%local_grid%size(z_index)+1-k, &
115  current_state%column_local_y, current_state%column_local_x) .gt. qlcrit) then
116  clbas_col=current_state%global_grid%configuration%vertical%zn(current_state%local_grid%size(z_index)+1-k)
117  end if
118 
119  !tql_loc(J,I)=tql_loc(J,I)+DZ(K)*RHON(K)*Q(J,K,IQL,current_state%liquid_water_mixing_ratio_index)
120  totqz=totqz+tempfac(k)*current_state%q(current_state%liquid_water_mixing_ratio_index)%data(k, &
121  current_state%column_local_y, current_state%column_local_x)*current_state%global_grid%configuration%vertical%zn(k)
122  end if
123 
124  if (current_state%water_vapour_mixing_ratio_index .gt. 0 .and. &
125  current_state%number_q_fields .ge. current_state%water_vapour_mixing_ratio_index) then
126  totqv=totqv+tempfac(k)*current_state%q(current_state%water_vapour_mixing_ratio_index)%data(k, &
127  current_state%column_local_y, current_state%column_local_x)
128  totql=totql+tempfac(k)*current_state%q(current_state%liquid_water_mixing_ratio_index)%data(k, &
129  current_state%column_local_y, current_state%column_local_x)
130  end if
131  end if
132  end do
133  if (cltop_col .gt. 0.0_default_precision) ncl_col=ncl_col+1
134  cltop_av=cltop_av+cltop_col
135  clbas_av=clbas_av+clbas_col
136 
137  end if
Here is the caller graph for this function:

Variable Documentation

◆ clbas

real(kind=default_precision) diagnostics_mod::clbas
private

Definition at line 18 of file diagnostics.F90.

◆ clbas_av

real(kind=default_precision) diagnostics_mod::clbas_av
private

Definition at line 18 of file diagnostics.F90.

◆ cltop

real(kind=default_precision) diagnostics_mod::cltop
private

Definition at line 18 of file diagnostics.F90.

◆ cltop_av

real(kind=default_precision) diagnostics_mod::cltop_av
private

Definition at line 18 of file diagnostics.F90.

◆ hqlmax

integer diagnostics_mod::hqlmax
private

Definition at line 16 of file diagnostics.F90.

◆ ncl_col

integer diagnostics_mod::ncl_col
private

Definition at line 16 of file diagnostics.F90.

◆ qlcrit

real(kind=default_precision) diagnostics_mod::qlcrit
private

Definition at line 18 of file diagnostics.F90.

◆ qlmax

real(kind=default_precision) diagnostics_mod::qlmax
private

Definition at line 18 of file diagnostics.F90.

◆ tempfac

real(kind=default_precision), dimension(:), allocatable diagnostics_mod::tempfac
private

Definition at line 17 of file diagnostics.F90.

17  real(kind=DEFAULT_PRECISION), dimension(:), allocatable :: tempfac, theta_tot

◆ theta_tot

real(kind=default_precision), dimension(:), allocatable diagnostics_mod::theta_tot
private

Definition at line 17 of file diagnostics.F90.

◆ total_points

integer diagnostics_mod::total_points
private

Definition at line 16 of file diagnostics.F90.

16  integer :: total_points, hqlmax, ncl_col

◆ totql

real(kind=default_precision) diagnostics_mod::totql
private

Definition at line 18 of file diagnostics.F90.

◆ totqv

real(kind=default_precision) diagnostics_mod::totqv
private

Definition at line 18 of file diagnostics.F90.

18  real(kind=DEFAULT_PRECISION) :: totqv, totql, totqz, wmax, wmin, qlmax, qlcrit, cltop_av, clbas_av, cltop, clbas

◆ totqz

real(kind=default_precision) diagnostics_mod::totqz
private

Definition at line 18 of file diagnostics.F90.

◆ wmax

real(kind=default_precision) diagnostics_mod::wmax
private

Definition at line 18 of file diagnostics.F90.

◆ wmin

real(kind=default_precision) diagnostics_mod::wmin
private

Definition at line 18 of file diagnostics.F90.