MONC
Functions/Subroutines | Variables
scalar_diagnostics_mod Module Reference

Functions/Subroutines

type(component_descriptor_type) function, public scalar_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 ncl_col
 
real(kind=default_precision), dimension(:), allocatable tempfac
 
real(kind=default_precision) qlcrit
 
real(kind=default_precision), dimension(:,:), allocatable vwp
 
real(kind=default_precision), dimension(:,:), allocatable lwp
 
real(kind=default_precision), dimension(:,:), allocatable wmax
 
real(kind=default_precision), dimension(:,:), allocatable wmin
 
real(kind=default_precision), dimension(:,:), allocatable qlmax
 
real(kind=default_precision), dimension(:,:), allocatable hqlmax
 
real(kind=default_precision), dimension(:,:), allocatable cltop
 
real(kind=default_precision), dimension(:,:), allocatable clbas
 
real(kind=default_precision), dimension(:,:), allocatable senhf
 
real(kind=default_precision), dimension(:,:), allocatable lathf
 

Function/Subroutine Documentation

◆ field_information_retrieval_callback()

subroutine scalar_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 195 of file scalar_diagnostics.F90.

195  type(model_state_type), target, intent(inout) :: current_state
196  character(len=*), intent(in) :: name
197  type(component_field_information_type), intent(out) :: field_information
198 
199  field_information%field_type=component_array_field_type
200  field_information%number_dimensions=2
201  field_information%dimension_sizes(1)=current_state%local_grid%size(y_index)
202  field_information%dimension_sizes(2)=current_state%local_grid%size(x_index)
203  field_information%data_type=component_double_data_type
204 
205  if (name .eq. "senhf_local") then
206  field_information%enabled=current_state%use_surface_boundary_conditions .and. current_state%th%active
207  else if (name .eq. "lathf_local") then
208  field_information%enabled=current_state%use_surface_boundary_conditions .and. &
209  current_state%water_vapour_mixing_ratio_index .gt. 0 .and. &
210  current_state%number_q_fields .ge. current_state%water_vapour_mixing_ratio_index
211  else if (name .eq. "cltop_local" .or. name .eq. "clbas_local") then
212  field_information%enabled=current_state%number_q_fields .gt. 0
213  else if (name .eq. "qlmax_local") then
214  field_information%enabled=current_state%number_q_fields .gt. 0 .and. current_state%liquid_water_mixing_ratio_index .gt. 0 &
215  .and. current_state%number_q_fields .ge. current_state%liquid_water_mixing_ratio_index
216  else if (name .eq. "vwp_local" .or. name .eq. "lwp_local") then
217  field_information%enabled=current_state%number_q_fields .gt. 0 .and. current_state%water_vapour_mixing_ratio_index .gt. 0 &
218  .and. current_state%number_q_fields .ge. current_state%water_vapour_mixing_ratio_index
219  else
220  field_information%enabled=.true.
221  end if
222 
Here is the caller graph for this function:

◆ field_value_retrieval_callback()

subroutine scalar_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 230 of file scalar_diagnostics.F90.

230  type(model_state_type), target, intent(inout) :: current_state
231  character(len=*), intent(in) :: name
232  type(component_field_value_type), intent(out) :: field_value
233 
234  integer :: i
235 
236  if (name .eq. "wmax_local") then
237  allocate(field_value%real_2d_array(current_state%local_grid%size(y_index), &
238  current_state%local_grid%size(x_index)))
239  field_value%real_2d_array(:,:)=wmax(:,:)
240  else if (name .eq. "wmin_local") then
241  allocate(field_value%real_2d_array(current_state%local_grid%size(y_index), &
242  current_state%local_grid%size(x_index)))
243  field_value%real_2d_array(:,:)=wmin(:,:)
244  else if (name .eq. "qlmax_local") then
245  allocate(field_value%real_2d_array(current_state%local_grid%size(y_index), &
246  current_state%local_grid%size(x_index)))
247  field_value%real_2d_array(:,:)=qlmax(:,:)
248  !else if (name .eq. "hqlmax_local") then
249  ! allocate(field_value%real_2d_array(current_state%local_grid%size(Y_INDEX), &
250  ! current_state%local_grid%size(X_INDEX)))
251  ! field_value%real_2d_array(:,:)=hqlmax(:,:)
252  else if (name .eq. "cltop_local") then
253  allocate(field_value%real_2d_array(current_state%local_grid%size(y_index), &
254  current_state%local_grid%size(x_index)))
255  field_value%real_2d_array(:,:)=cltop(:,:)
256  else if (name .eq. "clbas_local") then
257  allocate(field_value%real_2d_array(current_state%local_grid%size(y_index), &
258  current_state%local_grid%size(x_index)))
259  field_value%real_2d_array(:,:)=clbas(:,:)
260  else if (name .eq. "vwp_local") then
261  allocate(field_value%real_2d_array(current_state%local_grid%size(y_index), &
262  current_state%local_grid%size(x_index)))
263  field_value%real_2d_array(:,:)=vwp(:,:)
264  else if (name .eq. "lwp_local") then
265  allocate(field_value%real_2d_array(current_state%local_grid%size(y_index), &
266  current_state%local_grid%size(x_index)))
267  field_value%real_2d_array(:,:)=lwp(:,:)
268  else if (name .eq. "senhf_local") then
269  allocate(field_value%real_2d_array(current_state%local_grid%size(y_index), &
270  current_state%local_grid%size(x_index)))
271  field_value%real_2d_array(:,:)=senhf(:,:)
272  else if (name .eq. "lathf_local") then
273  allocate(field_value%real_2d_array(current_state%local_grid%size(y_index), &
274  current_state%local_grid%size(x_index)))
275  field_value%real_2d_array(:,:)=lathf(:,:)
276  end if
277 
Here is the caller graph for this function:

◆ initialisation_callback()

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

Definition at line 55 of file scalar_diagnostics.F90.

55  type(model_state_type), target, intent(inout) :: current_state
56 
57  integer :: k
58  integer :: y_size_local, x_size_local
59 
60  ! qlcrit declared in the config file
61  qlcrit=options_get_real(current_state%options_database, "qlcrit")
62 
63  y_size_local = current_state%local_grid%size(y_index)
64  x_size_local = current_state%local_grid%size(x_index)
65 
66  ! allocate scalar diagnostics as 2-D fields (horizontal slices) that
67  ! are published to the ioserver so that the ioserver can do the manipulation
68  ! to obtain the scalar field
69  allocate(vwp(y_size_local, x_size_local), lwp(y_size_local, x_size_local), &
70  wmax(y_size_local, x_size_local), wmin(y_size_local, x_size_local), &
71  qlmax(y_size_local, x_size_local), hqlmax(y_size_local, x_size_local), &
72  cltop(y_size_local, x_size_local), clbas(y_size_local, x_size_local), &
73  senhf(y_size_local, x_size_local),lathf(y_size_local, x_size_local))
74 
75  allocate(tempfac(current_state%local_grid%size(z_index)))
76 
77  do k=2, current_state%local_grid%size(z_index)
78  tempfac(k)=current_state%global_grid%configuration%vertical%dz(k)*&
79  current_state%global_grid%configuration%vertical%rhon(k)
80  end do
81 
Here is the caller graph for this function:

◆ scalar_diagnostics_get_descriptor()

type(component_descriptor_type) function, public scalar_diagnostics_mod::scalar_diagnostics_get_descriptor ( )

Provides the component descriptor for the core to register.

Returns
The descriptor describing this component

Definition at line 30 of file scalar_diagnostics.F90.

30  scalar_diagnostics_get_descriptor%name="scalar_diagnostics"
31  scalar_diagnostics_get_descriptor%version=0.1
32 
33  scalar_diagnostics_get_descriptor%initialisation=>initialisation_callback
34  scalar_diagnostics_get_descriptor%timestep=>timestep_callback
35 
36  scalar_diagnostics_get_descriptor%field_value_retrieval=>field_value_retrieval_callback
37  scalar_diagnostics_get_descriptor%field_information_retrieval=>field_information_retrieval_callback
38  allocate(scalar_diagnostics_get_descriptor%published_fields(10))
39 
40  scalar_diagnostics_get_descriptor%published_fields(1)="vwp_local"
41  scalar_diagnostics_get_descriptor%published_fields(2)="lwp_local"
42  scalar_diagnostics_get_descriptor%published_fields(3)="qlmax_local"
43  scalar_diagnostics_get_descriptor%published_fields(4)="hqlmax_local"
44  scalar_diagnostics_get_descriptor%published_fields(5)="cltop_local"
45  scalar_diagnostics_get_descriptor%published_fields(6)="clbas_local"
46  scalar_diagnostics_get_descriptor%published_fields(7)="wmax_local"
47  scalar_diagnostics_get_descriptor%published_fields(8)="wmin_local"
48  scalar_diagnostics_get_descriptor%published_fields(9)="senhf_local"
49  scalar_diagnostics_get_descriptor%published_fields(10)="lathf_local"
50 
51 
Here is the call graph for this function:

◆ timestep_callback()

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

Definition at line 85 of file scalar_diagnostics.F90.

85  type(model_state_type), target, intent(inout) :: current_state
86 
87  integer :: k, i
88  integer :: current_y_index, current_x_index, target_x_index, target_y_index
89 
90  current_y_index=current_state%column_local_y
91  current_x_index=current_state%column_local_x
92  target_y_index=current_y_index-current_state%local_grid%halo_size(y_index)
93  target_x_index=current_x_index-current_state%local_grid%halo_size(x_index)
94 
95  if (current_state%first_timestep_column) then
96  ! water vapour path for each column
97  vwp(:,:)=0.0
98  ! liquid water path for each column
99  lwp(:,:)=0.0
100  ! maximum vertical velocity for each column
101  wmax(:,:)=0.0
102  ! minimum vertical velocity for each column
103  wmin(:,:)=0.0
104  ! maximum liquid water content in a column
105  qlmax(:,:)=0.0
106  ! the height of the maximum liquid water content in a column
107  hqlmax(:,:)=0.0
108  ! cloud top height where liqud water content is greater than qlcrit
109  cltop(:,:)=0.0
110  ! minimum cloud base where liquid water content is greater than q;crit
111  clbas(:,:)=0.0
112  ! surface sensible heat flux
113  senhf(:,:)=0.0
114  ! surface latent heat flux
115  lathf(:,:)=0.0
116  end if
117 
118  if (.not. current_state%halo_column) then
119  if (current_state%number_q_fields .gt. 0) then
120  !
121  ! calculate the lwc maximum and height of lwc max for each column
122  !
123  if (current_state%liquid_water_mixing_ratio_index .gt. 0 .and. &
124  current_state%number_q_fields .ge. current_state%liquid_water_mixing_ratio_index) then
125  qlmax(target_y_index, target_x_index) = &
126  maxval(current_state%q(current_state%liquid_water_mixing_ratio_index)%data &
127  (:,current_y_index, current_x_index))
128  !hqlmax(current_y_index, current_x_index) = &
129  ! current_state%global_grid%configuration%vertical% &
130  ! zn(maxloc(current_state%q(current_state%liquid_water_mixing_ratio_index)%data &
131  ! (:,current_y_index, current_x_index)))
132  end if
133  !
134  ! calculate the cloud top maximum and minimum for each column
135  !
136  do k = 2, current_state%local_grid%size(z_index)
137  if (current_state%q(current_state%liquid_water_mixing_ratio_index)%data(k, &
138  current_y_index, current_x_index) .gt. qlcrit) then
139  cltop(target_y_index, target_x_index) = &
140  current_state%global_grid%configuration%vertical%zn(k)
141  endif
142  if (current_state%q(current_state%liquid_water_mixing_ratio_index)%data( &
143  current_state%local_grid%size(z_index)+1-k, current_y_index, current_x_index) .gt. &
144  qlcrit) then
145  clbas(target_y_index, target_x_index)= &
146  current_state%global_grid%configuration%vertical%zn(current_state%local_grid%size(z_index)+1-k)
147  end if
148  enddo
149  !
150  ! calculate the vapour and liquid water path
151  !
152  if (current_state%water_vapour_mixing_ratio_index .gt. 0 .and. &
153  current_state%number_q_fields .ge. current_state%water_vapour_mixing_ratio_index) then
154  do k = 2, current_state%local_grid%size(z_index)
155  vwp(target_y_index, target_x_index)=vwp(target_y_index, target_x_index) &
156  +tempfac(k)*current_state%q(current_state%water_vapour_mixing_ratio_index)%data(k, &
157  current_y_index, current_x_index)
158  lwp(target_y_index, target_x_index)=lwp(target_y_index, target_x_index) &
159  +tempfac(k)*current_state%q(current_state%liquid_water_mixing_ratio_index)%data(k, &
160  current_y_index, current_x_index)
161  enddo
162  end if
163  end if
164 
165  ! surface flux diagnostics
166  if (current_state%use_surface_boundary_conditions) then
167  if (current_state%water_vapour_mixing_ratio_index .gt. 0 .and. &
168  current_state%number_q_fields .ge. current_state%water_vapour_mixing_ratio_index) then
169  lathf(target_y_index, target_x_index)= (current_state%diff_coefficient%data(1, current_y_index, current_x_index) * &
170  current_state%global_grid%configuration%vertical%rdzn(2) * &
171  (current_state%q(current_state%water_vapour_mixing_ratio_index)%data(1,current_y_index,current_x_index) - &
172  current_state%q(current_state%water_vapour_mixing_ratio_index)%data(2,current_y_index,current_x_index))) &
173  * rlvap * current_state%global_grid%configuration%vertical%rhon(1)
174  endif
175 
176  if (current_state%th%active) then
177  senhf(target_y_index, target_x_index)=(current_state%diff_coefficient%data(1, current_y_index, current_x_index) &
178  * current_state%global_grid%configuration%vertical%rdzn(2) &
179  * (current_state%th%data(1, current_y_index, current_x_index) &
180  - current_state%th%data(2, current_y_index, current_x_index) &
181  - current_state%global_grid%configuration%vertical%dthref(1))) &
182  * current_state%global_grid%configuration%vertical%rhon(1)*cp
183  endif
184  endif
185  wmax(target_y_index, target_x_index)=maxval(current_state%w%data(:, current_y_index, current_x_index))
186  wmin(target_y_index, target_x_index)=minval(current_state%w%data(:, current_y_index, current_x_index))
187  end if
Here is the caller graph for this function:

Variable Documentation

◆ clbas

real(kind=default_precision), dimension(:,:), allocatable scalar_diagnostics_mod::clbas
private

Definition at line 20 of file scalar_diagnostics.F90.

◆ cltop

real(kind=default_precision), dimension(:,:), allocatable scalar_diagnostics_mod::cltop
private

Definition at line 20 of file scalar_diagnostics.F90.

◆ hqlmax

real(kind=default_precision), dimension(:,:), allocatable scalar_diagnostics_mod::hqlmax
private

Definition at line 20 of file scalar_diagnostics.F90.

◆ lathf

real(kind=default_precision), dimension(:,:), allocatable scalar_diagnostics_mod::lathf
private

Definition at line 20 of file scalar_diagnostics.F90.

◆ lwp

real(kind=default_precision), dimension(:,:), allocatable scalar_diagnostics_mod::lwp
private

Definition at line 20 of file scalar_diagnostics.F90.

◆ ncl_col

integer scalar_diagnostics_mod::ncl_col
private

Definition at line 17 of file scalar_diagnostics.F90.

◆ qlcrit

real(kind=default_precision) scalar_diagnostics_mod::qlcrit
private

Definition at line 19 of file scalar_diagnostics.F90.

19  real(kind=DEFAULT_PRECISION) :: qlcrit

◆ qlmax

real(kind=default_precision), dimension(:,:), allocatable scalar_diagnostics_mod::qlmax
private

Definition at line 20 of file scalar_diagnostics.F90.

◆ senhf

real(kind=default_precision), dimension(:,:), allocatable scalar_diagnostics_mod::senhf
private

Definition at line 20 of file scalar_diagnostics.F90.

◆ tempfac

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

Definition at line 18 of file scalar_diagnostics.F90.

18  real(kind=DEFAULT_PRECISION), dimension(:), allocatable :: tempfac

◆ total_points

integer scalar_diagnostics_mod::total_points
private

Definition at line 17 of file scalar_diagnostics.F90.

17  integer :: total_points, ncl_col

◆ vwp

real(kind=default_precision), dimension(:,:), allocatable scalar_diagnostics_mod::vwp
private

Definition at line 20 of file scalar_diagnostics.F90.

20  real(kind=DEFAULT_PRECISION), dimension(:,:), allocatable :: vwp, lwp, wmax, wmin, &
21  qlmax, hqlmax, cltop, clbas, senhf, lathf

◆ wmax

real(kind=default_precision), dimension(:,:), allocatable scalar_diagnostics_mod::wmax
private

Definition at line 20 of file scalar_diagnostics.F90.

◆ wmin

real(kind=default_precision), dimension(:,:), allocatable scalar_diagnostics_mod::wmin
private

Definition at line 20 of file scalar_diagnostics.F90.