MONC
profile_diagnostics.F90
Go to the documentation of this file.
5  use grids_mod, only : z_index, y_index, x_index
7  use state_mod, only : model_state_type
10  use saturation_mod, only: qsaturation
11 
12  implicit none
13 
14 #ifndef TEST_MODE
15  private
16 #endif
17 
18  integer :: total_points, iqv, iql
19  real(kind=DEFAULT_PRECISION), dimension(:), allocatable :: &
24  real(kind=DEFAULT_PRECISION) :: qlcrit
25 
27 
28 contains
29 
33  profile_diagnostics_get_descriptor%name="profile_diagnostics"
35 
38 
41  allocate(profile_diagnostics_get_descriptor%published_fields(2*22))
42 
43  profile_diagnostics_get_descriptor%published_fields(1)="theta_total_local"
44  profile_diagnostics_get_descriptor%published_fields(2)="vapour_mmr_total_local"
45  profile_diagnostics_get_descriptor%published_fields(3)="liquid_mmr_total_local"
46  profile_diagnostics_get_descriptor%published_fields(4)="u_wind_total_local"
47  profile_diagnostics_get_descriptor%published_fields(5)="uu_total_local"
48  profile_diagnostics_get_descriptor%published_fields(6)="v_wind_total_local"
49  profile_diagnostics_get_descriptor%published_fields(7)="vv_total_local"
50  profile_diagnostics_get_descriptor%published_fields(8)="ww_total_local"
51  profile_diagnostics_get_descriptor%published_fields(9)="w_wind_total_local"
52  profile_diagnostics_get_descriptor%published_fields(10)="rh_total_local"
53  profile_diagnostics_get_descriptor%published_fields(11)="thref_local"
54  profile_diagnostics_get_descriptor%published_fields(12)="prefn_local"
55  profile_diagnostics_get_descriptor%published_fields(13)="rho_local"
56  profile_diagnostics_get_descriptor%published_fields(14)="rhon_local"
57  profile_diagnostics_get_descriptor%published_fields(15)="thinit_local"
58  profile_diagnostics_get_descriptor%published_fields(16)="uw_total_local"
59  profile_diagnostics_get_descriptor%published_fields(17)="vw_total_local"
60  profile_diagnostics_get_descriptor%published_fields(18)="wtheta_total_local"
61  profile_diagnostics_get_descriptor%published_fields(19)="th2_total_local"
62  profile_diagnostics_get_descriptor%published_fields(20)="wq_total_local"
63 
64  profile_diagnostics_get_descriptor%published_fields(21)="uinit_local"
65  profile_diagnostics_get_descriptor%published_fields(22)="vinit_local"
66 
67 ! =====================================================
68 ! 2nd, provisionally instantaneous, stream
69  profile_diagnostics_get_descriptor%published_fields(22+1)="i_theta_total_local"
70  profile_diagnostics_get_descriptor%published_fields(22+2)="i_vapour_mmr_total_local"
71  profile_diagnostics_get_descriptor%published_fields(22+3)="i_liquid_mmr_total_local"
72  profile_diagnostics_get_descriptor%published_fields(22+4)="i_u_wind_total_local"
73  profile_diagnostics_get_descriptor%published_fields(22+4)="i_u_wind_total_local"
74  profile_diagnostics_get_descriptor%published_fields(22+5)="i_uu_total_local"
75  profile_diagnostics_get_descriptor%published_fields(22+6)="i_v_wind_total_local"
76  profile_diagnostics_get_descriptor%published_fields(22+7)="i_vv_total_local"
77  profile_diagnostics_get_descriptor%published_fields(22+8)="i_ww_total_local"
78  profile_diagnostics_get_descriptor%published_fields(22+9)="i_w_wind_total_local"
79  profile_diagnostics_get_descriptor%published_fields(22+10)="i_rh_total_local"
80  profile_diagnostics_get_descriptor%published_fields(22+11)="i_thref_local"
81  profile_diagnostics_get_descriptor%published_fields(22+12)="i_prefn_local"
82  profile_diagnostics_get_descriptor%published_fields(22+13)="i_rho_local"
83  profile_diagnostics_get_descriptor%published_fields(22+14)="i_rhon_local"
84  profile_diagnostics_get_descriptor%published_fields(22+15)="i_thinit_local"
85 
86  profile_diagnostics_get_descriptor%published_fields(22+16)="i_uw_total_local"
87  profile_diagnostics_get_descriptor%published_fields(22+17)="i_vw_total_local"
88  profile_diagnostics_get_descriptor%published_fields(22+18)="i_wtheta_total_local"
89  profile_diagnostics_get_descriptor%published_fields(22+19)="i_th2_total_local"
90  profile_diagnostics_get_descriptor%published_fields(22+20)="i_wq_total_local"
91 
92  profile_diagnostics_get_descriptor%published_fields(22+21)="i_uinit_local"
93  profile_diagnostics_get_descriptor%published_fields(22+22)="i_vinit_local"
94 
96 
97  subroutine initialisation_callback(current_state)
98  type(model_state_type), target, intent(inout) :: current_state
99 
100  integer :: k
101 
102  total_points=(current_state%local_grid%size(y_index) * current_state%local_grid%size(x_index))
103 
104  ! allocate local arrays for the horizontal wind averages
105  allocate(u_wind_tot(current_state%local_grid%size(z_index)) &
106  , v_wind_tot(current_state%local_grid%size(z_index)) &
107  , ww_tot(current_state%local_grid%size(z_index)) &
108  , w_wind_tot(current_state%local_grid%size(z_index)) &
109  , prefn(current_state%local_grid%size(z_index)) &
110  , rho(current_state%local_grid%size(z_index)) &
111  , rhon(current_state%local_grid%size(z_index)) &
112  , uinit(current_state%local_grid%size(z_index)) &
113  , vinit(current_state%local_grid%size(z_index)))
114  allocate(wq_tot(current_state%local_grid%size(z_index)) &
115  , wtheta_tot(current_state%local_grid%size(z_index)) &
116  , uw_tot(current_state%local_grid%size(z_index)) &
117  , vw_tot(current_state%local_grid%size(z_index)) &
118  , th2_tot(current_state%local_grid%size(z_index)) )
119 
120  if (allocated(current_state%global_grid%configuration%vertical%olubar)) then
121  allocate(uprime_tot(current_state%local_grid%size(z_index)))
122  end if
123 
124  if (allocated(current_state%global_grid%configuration%vertical%olvbar)) then
125  allocate(vprime_tot(current_state%local_grid%size(z_index)))
126  end if
127 
128  if (current_state%th%active) then
129  allocate(theta_tot(current_state%local_grid%size(z_index)) &
130  , thref(current_state%local_grid%size(z_index)) &
131  , thinit(current_state%local_grid%size(z_index)))
132  endif
133 
134  if (.not. current_state%passive_q .and. current_state%number_q_fields .gt. 0) then
135  iqv=get_q_index(standard_q_names%VAPOUR, 'profile_diags')
136  iql=get_q_index(standard_q_names%CLOUD_LIQUID_MASS, 'profile_diags')
137  qlcrit=options_get_real(current_state%options_database, "qlcrit")
138  allocate(qv_tot(current_state%local_grid%size(z_index)) &
139  , ql_tot(current_state%local_grid%size(z_index)) )
140  if (current_state%th%active) &
141  allocate(rh_tot(current_state%local_grid%size(z_index)))
142  endif
143 
144  end subroutine initialisation_callback
145 
146  subroutine timestep_callback(current_state)
147  type(model_state_type), target, intent(inout) :: current_state
148 
149  integer :: k, i
150  real(kind=DEFAULT_PRECISION) :: cltop_col, clbas_col, qv, qc, TdegK, Pmb &
151  , qs, exner
152  real(kind=DEFAULT_PRECISION) :: uprime_w_local, vprime_w_local &
153  , thprime_w_local, qprime_w_local
154 
155  if (current_state%first_timestep_column) then
156  u_wind_tot(:) = 0.0_default_precision
157  if (allocated(uprime_tot)) uprime_tot(:) = 0.0_default_precision
158  v_wind_tot(:) = 0.0_default_precision
159  if (allocated(vprime_tot)) vprime_tot(:) = 0.0_default_precision
160  w_wind_tot(:) = 0.0_default_precision
161  ww_tot(:) = 0.0_default_precision
162 
163  wq_tot(:) = 0.0_default_precision
164  wtheta_tot(:) = 0.0_default_precision
165  uw_tot(:) = 0.0_default_precision
166  vw_tot(:) = 0.0_default_precision
167  th2_tot(:) = 0.0_default_precision
168 
169  if (current_state%th%active) then
170  theta_tot(:)=0.0_default_precision
171  endif
172  if (.not. current_state%passive_q .and. &
173  current_state%number_q_fields .gt. 0) then
174  qv_tot(:)=0.0_default_precision
175  ql_tot(:)=0.0_default_precision
176  if (current_state%th%active) &
177  rh_tot(:) = 0.0_default_precision
178  endif
179  end if
180  if (.not. current_state%halo_column) then
181  ! work out the sum of u and v wind over local domaini
182  do k=1, current_state%local_grid%size(z_index)
183  u_wind_tot(k) = u_wind_tot(k) + &
184  (current_state%u%data(k,current_state%column_local_y,current_state%column_local_x) &
185  + current_state%ugal)
186  if (allocated(uprime_tot)) then
187  uprime_tot(k) = uprime_tot(k) + &
188  ((current_state%u%data(k,current_state%column_local_y,current_state%column_local_x) &
189  - (current_state%global_grid%configuration%vertical%olubar(k) - current_state%ugal))**2.)
190  end if
191  v_wind_tot(k) = v_wind_tot(k) + &
192  (current_state%v%data(k,current_state%column_local_y,current_state%column_local_x) &
193  + current_state%vgal)
194  if (allocated(vprime_tot)) then
195  vprime_tot(k) = vprime_tot(k) + &
196  ((current_state%v%data(k,current_state%column_local_y,current_state%column_local_x) &
197  - (current_state%global_grid%configuration%vertical%olvbar(k) - current_state%vgal))**2.)
198  end if
199  ww_tot(k) = ww_tot(k) + &
200  (current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)**2.)
201  w_wind_tot(k) = w_wind_tot(k) + &
202  (current_state%w%data(k,current_state%column_local_y,current_state%column_local_x))
203  enddo
204 
205 ! <u'w'> and <v'w'> are on w-points, so we interpolate u and v both horizontally and vertically.
206  do k=1, current_state%local_grid%size(z_index)-1
207  uprime_w_local = &
208  0.25 * ( current_state%u%data(k,current_state%column_local_y,current_state%column_local_x) + &
209  current_state%u%data(k,current_state%column_local_y,current_state%column_local_x-1) + &
210  current_state%u%data(k+1,current_state%column_local_y,current_state%column_local_x) + &
211  current_state%u%data(k+1,current_state%column_local_y,current_state%column_local_x-1) ) + &
212  current_state%ugal
213  if (allocated(current_state%global_grid%configuration%vertical%olubar)) &
214  uprime_w_local = uprime_w_local - &
215  0.5 * ( current_state%global_grid%configuration%vertical%olubar(k) + &
216  current_state%global_grid%configuration%vertical%olubar(k+1) )
217  vprime_w_local = &
218  0.25 * ( current_state%v%data(k,current_state%column_local_y,current_state%column_local_x) + &
219  current_state%v%data(k,current_state%column_local_y-1,current_state%column_local_x) + &
220  current_state%v%data(k+1,current_state%column_local_y,current_state%column_local_x) + &
221  current_state%v%data(k+1,current_state%column_local_y-1,current_state%column_local_x) ) + &
222  current_state%vgal
223  if (allocated(current_state%global_grid%configuration%vertical%olvbar)) &
224  vprime_w_local = vprime_w_local - &
225  0.5 * ( current_state%global_grid%configuration%vertical%olvbar(k) + &
226  current_state%global_grid%configuration%vertical%olvbar(k+1) )
227  uw_tot(k) = uw_tot(k) + uprime_w_local * &
228  current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)
229  vw_tot(k) = vw_tot(k) + vprime_w_local * &
230  current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)
231  enddo
232 
233  if (current_state%th%active) then
234  do k=1, current_state%local_grid%size(z_index)
235  theta_tot(k) = theta_tot(k) + &
236  (current_state%th%data(k,current_state%column_local_y,current_state%column_local_x) &
237  + current_state%global_grid%configuration%vertical%thref(k))
238  th2_tot(k) = th2_tot(k) + &
239  (current_state%th%data(k,current_state%column_local_y,current_state%column_local_x) - &
240  current_state%global_grid%configuration%vertical%olthbar(k) )**2
241  enddo
242 ! <w'theta'> is on w-levels, so theta is interpolated to w-levels.
243  do k=1, current_state%local_grid%size(z_index)-1
244  thprime_w_local = 0.5 * ( &
245  (current_state%th%data(k,current_state%column_local_y,current_state%column_local_x) - &
246  current_state%global_grid%configuration%vertical%olthbar(k)) + &
247  (current_state%th%data(k+1,current_state%column_local_y,current_state%column_local_x) - &
248  current_state%global_grid%configuration%vertical%olthbar(k+1)) )
249  wtheta_tot(k) = wtheta_tot(k) + &
250  current_state%w%data(k,current_state%column_local_y,current_state%column_local_x) * &
251  thprime_w_local
252  enddo
253  endif
254  if (current_state%th%active .and. .not. current_state%passive_q .and. current_state%number_q_fields .gt. 0) then
255  do k=1, current_state%local_grid%size(z_index)
256  qv_tot(k) = qv_tot(k) + (current_state%q(iqv)%data(k,current_state%column_local_y,current_state%column_local_x))
257  ql_tot(k) = ql_tot(k) + (current_state%q(iql)%data(k,current_state%column_local_y,current_state%column_local_x))
258  ! temporary code for RH calculation
259  exner = current_state%global_grid%configuration%vertical%rprefrcp(k)
260  pmb = (current_state%global_grid%configuration%vertical%prefn(k)/100.)
261  qv = current_state%q(iqv)%data(k, current_state%column_local_y,current_state%column_local_x)
262  qc = current_state%q(iql)%data(k, current_state%column_local_y,current_state%column_local_x)
263  tdegk = (current_state%th%data(k,current_state%column_local_y,current_state%column_local_x) &
264  + current_state%global_grid%configuration%vertical%thref(k))*exner
265  qs = qsaturation(tdegk, pmb)
266  rh_tot(k) = rh_tot(k) + (qv/qs)
267  enddo
268  do k=1, current_state%local_grid%size(z_index)-1
269  qprime_w_local = 0.5 * ( &
270  (current_state%q(current_state%liquid_water_mixing_ratio_index)% &
271  data(k,current_state%column_local_y,current_state%column_local_x) - &
272  current_state%global_grid%configuration%vertical% &
273  olqbar(k,current_state%liquid_water_mixing_ratio_index)) + &
274  (current_state%q(current_state%liquid_water_mixing_ratio_index)% &
275  data(k+1,current_state%column_local_y,current_state%column_local_x) - &
276  current_state%global_grid%configuration%vertical% &
277  olqbar(k+1,current_state%liquid_water_mixing_ratio_index)) )
278  wq_tot(k) = wq_tot(k) + qprime_w_local * &
279  current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)
280  end do
281  endif
282  endif
283  end subroutine timestep_callback
284 
289  subroutine field_information_retrieval_callback(current_state, name, field_information)
290  type(model_state_type), target, intent(inout) :: current_state
291  character(len=*), intent(in) :: name
292  type(component_field_information_type), intent(out) :: field_information
293 
294  field_information%field_type=component_array_field_type
295  field_information%number_dimensions=1
296  field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
297  field_information%data_type=component_double_data_type
298  if (name .eq. "theta_total_local") then
299  field_information%enabled=current_state%th%active
300  else if (name .eq. "vapour_mmr_total_local" .or. name .eq. "liquid_mmr_total_local") then
301  field_information%enabled=.not. current_state%passive_q .and. current_state%number_q_fields .gt. 0
302  else if (name .eq. "rh_total_local") then
303  field_information%enabled=current_state%th%active .and. .not. current_state%passive_q .and. &
304  current_state%number_q_fields .gt. 0
305  else if (name .eq. "uu_total_local") then
306  field_information%enabled=allocated(uprime_tot)
307  else if (name .eq. "vv_total_local") then
308  field_information%enabled=allocated(vprime_tot)
309  else if (name .eq. "wtheta_total_local") then
310  field_information%enabled=current_state%th%active
311  else if (name .eq. "wq_total_local") then
312  field_information%enabled= (.not.current_state%passive_q) .and. &
313  (current_state%liquid_water_mixing_ratio_index > 0)
314  else if (name .eq. "th2_total_local") then
315  field_information%enabled=current_state%th%active
316 ! ========================================================================
317 ! 2nd stream
318  else if (name .eq. "i_theta_total_local") then
319  field_information%enabled=current_state%th%active
320  else if (name .eq. "i_vapour_mmr_total_local" .or. name .eq. "i_liquid_mmr_total_local") then
321  field_information%enabled=.not. current_state%passive_q .and. current_state%number_q_fields .gt. 0
322  else if (name .eq. "i_rh_total_local") then
323  field_information%enabled=current_state%th%active .and. .not. current_state%passive_q .and. &
324  current_state%number_q_fields .gt. 0
325  else if (name .eq. "i_uu_total_local") then
326  field_information%enabled=allocated(uprime_tot)
327  else if (name .eq. "i_vv_total_local") then
328  field_information%enabled=allocated(vprime_tot)
329  else if (name .eq. "i_wtheta_total_local") then
330  field_information%enabled=current_state%th%active
331  else if (name .eq. "i_wq_total_local") then
332  field_information%enabled= (.not.current_state%passive_q) .and. &
333  (current_state%liquid_water_mixing_ratio_index > 0)
334  else if (name .eq. "i_th2_total_local") then
335  field_information%enabled=current_state%th%active
336 ! ========================================================================
337  else
338  field_information%enabled=.true.
339  end if
341 
346  subroutine field_value_retrieval_callback(current_state, name, field_value)
347  type(model_state_type), target, intent(inout) :: current_state
348  character(len=*), intent(in) :: name
349  type(component_field_value_type), intent(out) :: field_value
350 
351  integer :: k
352 
353  if (name .eq. "prefn_local") then
354  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
355  do k = 1, current_state%local_grid%size(z_index)
356  field_value%real_1d_array(k)=current_state%global_grid%configuration%vertical%prefn(k)
357  enddo
358  else if (name .eq. "rho_local") then
359  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
360  do k = 1, current_state%local_grid%size(z_index)
361  field_value%real_1d_array(k)=current_state%global_grid%configuration%vertical%rho(k)
362  enddo
363  else if (name .eq. "rhon_local") then
364  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
365  do k = 1, current_state%local_grid%size(z_index)
366  field_value%real_1d_array(k)=current_state%global_grid%configuration%vertical%rhon(k)
367  enddo
368  else if (name .eq. "thref_local") then
369  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
370  do k = 1, current_state%local_grid%size(z_index)
371  field_value%real_1d_array(k)=current_state%global_grid%configuration%vertical%thref(k)
372  enddo
373  else if (name .eq. "thinit_local") then
374  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
375  do k = 1, current_state%local_grid%size(z_index)
376  field_value%real_1d_array(k)=current_state%global_grid%configuration%vertical%theta_init(k)
377  enddo
378  elseif (name .eq. "u_wind_total_local") then
379  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
380  do k = 1, current_state%local_grid%size(z_index)
381  field_value%real_1d_array(k)=u_wind_tot(k)
382  enddo
383  else if (name .eq. "uu_total_local") then
384  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
385  do k = 1, current_state%local_grid%size(z_index)
386  field_value%real_1d_array(k)=uprime_tot(k)
387  enddo
388  else if (name .eq. "v_wind_total_local") then
389  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
390  do k = 1, current_state%local_grid%size(z_index)
391  field_value%real_1d_array(k)=v_wind_tot(k)
392  enddo
393  else if (name .eq. "vv_total_local") then
394  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
395  do k = 1, current_state%local_grid%size(z_index)
396  field_value%real_1d_array(k)=vprime_tot(k)
397  enddo
398  else if (name .eq. "ww_total_local") then
399  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
400  do k = 1, current_state%local_grid%size(z_index)
401  field_value%real_1d_array(k)=ww_tot(k)
402  enddo
403  else if (name .eq. "theta_total_local") then
404  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
405  do k = 1, current_state%local_grid%size(z_index)
406  field_value%real_1d_array(k)=theta_tot(k)
407  enddo
408  else if (name .eq. "vapour_mmr_total_local") then
409  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
410  do k = 1, current_state%local_grid%size(z_index)
411  field_value%real_1d_array(k)=qv_tot(k)
412  enddo
413  else if (name .eq. "liquid_mmr_total_local") then
414  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
415  do k = 1, current_state%local_grid%size(z_index)
416  field_value%real_1d_array(k)=ql_tot(k)
417  enddo
418  else if (name .eq. "w_wind_total_local") then
419  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
420  do k = 1, current_state%local_grid%size(z_index)
421  field_value%real_1d_array(k)=w_wind_tot(k)
422  enddo
423  else if (name .eq. "rh_total_local") then
424  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
425  do k = 1, current_state%local_grid%size(z_index)
426  field_value%real_1d_array(k)=rh_tot(k)
427  enddo
428  else if (name .eq. "wq_total_local") then
429  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
430  do k = 1, current_state%local_grid%size(z_index)
431  field_value%real_1d_array(k)=wq_tot(k)
432  enddo
433  else if (name .eq. "wtheta_total_local") then
434  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
435  do k = 1, current_state%local_grid%size(z_index)
436  field_value%real_1d_array(k)=wtheta_tot(k)
437  enddo
438  else if (name .eq. "uw_total_local") then
439  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
440  do k = 1, current_state%local_grid%size(z_index)
441  field_value%real_1d_array(k)=uw_tot(k)
442  enddo
443  else if (name .eq. "vw_total_local") then
444  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
445  do k = 1, current_state%local_grid%size(z_index)
446  field_value%real_1d_array(k)=vw_tot(k)
447  enddo
448  else if (name .eq. "th2_total_local") then
449  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
450  do k = 1, current_state%local_grid%size(z_index)
451  field_value%real_1d_array(k)=th2_tot(k)
452  enddo
453 ! =====================================================
454 ! 2nd stream
455  else if (name .eq. "i_prefn_local") then
456  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
457  do k = 1, current_state%local_grid%size(z_index)
458  field_value%real_1d_array(k)=current_state%global_grid%configuration%vertical%prefn(k)
459  enddo
460  else if (name .eq. "i_rho_local") then
461  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
462  do k = 1, current_state%local_grid%size(z_index)
463  field_value%real_1d_array(k)=current_state%global_grid%configuration%vertical%rho(k)
464  enddo
465  else if (name .eq. "i_rhon_local") then
466  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
467  do k = 1, current_state%local_grid%size(z_index)
468  field_value%real_1d_array(k)=current_state%global_grid%configuration%vertical%rhon(k)
469  enddo
470  else if (name .eq. "i_thref_local") then
471  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
472  do k = 1, current_state%local_grid%size(z_index)
473  field_value%real_1d_array(k)=current_state%global_grid%configuration%vertical%thref(k)
474  enddo
475  else if (name .eq. "i_thinit_local") then
476  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
477  do k = 1, current_state%local_grid%size(z_index)
478  field_value%real_1d_array(k)=current_state%global_grid%configuration%vertical%theta_init(k)
479  enddo
480  else if (name .eq. "i_uinit_local") then
481  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
482  do k = 1, current_state%local_grid%size(z_index)
483  field_value%real_1d_array(k)=current_state%global_grid%configuration%vertical%u_init(k)
484  enddo
485  else if (name .eq. "i_vinit_local") then
486  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
487  do k = 1, current_state%local_grid%size(z_index)
488  field_value%real_1d_array(k)=current_state%global_grid%configuration%vertical%v_init(k)
489  enddo
490  elseif (name .eq. "i_u_wind_total_local") then
491  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
492  do k = 1, current_state%local_grid%size(z_index)
493  field_value%real_1d_array(k)=u_wind_tot(k)
494  enddo
495  else if (name .eq. "i_uu_total_local") then
496  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
497  do k = 1, current_state%local_grid%size(z_index)
498  field_value%real_1d_array(k)=uprime_tot(k)
499  enddo
500  else if (name .eq. "i_v_wind_total_local") then
501  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
502  do k = 1, current_state%local_grid%size(z_index)
503  field_value%real_1d_array(k)=v_wind_tot(k)
504  enddo
505  else if (name .eq. "i_vv_total_local") then
506  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
507  do k = 1, current_state%local_grid%size(z_index)
508  field_value%real_1d_array(k)=vprime_tot(k)
509  enddo
510  else if (name .eq. "i_ww_total_local") then
511  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
512  do k = 1, current_state%local_grid%size(z_index)
513  field_value%real_1d_array(k)=ww_tot(k)
514  enddo
515  else if (name .eq. "i_theta_total_local") then
516  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
517  do k = 1, current_state%local_grid%size(z_index)
518  field_value%real_1d_array(k)=theta_tot(k)
519  enddo
520  else if (name .eq. "i_vapour_mmr_total_local") then
521  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
522  do k = 1, current_state%local_grid%size(z_index)
523  field_value%real_1d_array(k)=qv_tot(k)
524  enddo
525  else if (name .eq. "i_liquid_mmr_total_local") then
526  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
527  do k = 1, current_state%local_grid%size(z_index)
528  field_value%real_1d_array(k)=ql_tot(k)
529  enddo
530  else if (name .eq. "i_w_wind_total_local") then
531  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
532  do k = 1, current_state%local_grid%size(z_index)
533  field_value%real_1d_array(k)=w_wind_tot(k)
534  enddo
535  else if (name .eq. "i_rh_total_local") then
536  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
537  do k = 1, current_state%local_grid%size(z_index)
538  field_value%real_1d_array(k)=rh_tot(k)
539  enddo
540  else if (name .eq. "i_wq_total_local") then
541  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
542  do k = 1, current_state%local_grid%size(z_index)
543  field_value%real_1d_array(k)=wq_tot(k)
544  enddo
545  else if (name .eq. "i_wtheta_total_local") then
546  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
547  do k = 1, current_state%local_grid%size(z_index)
548  field_value%real_1d_array(k)=wtheta_tot(k)
549  enddo
550  else if (name .eq. "i_uw_total_local") then
551  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
552  do k = 1, current_state%local_grid%size(z_index)
553  field_value%real_1d_array(k)=uw_tot(k)
554  enddo
555  else if (name .eq. "i_vw_total_local") then
556  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
557  do k = 1, current_state%local_grid%size(z_index)
558  field_value%real_1d_array(k)=vw_tot(k)
559  enddo
560  else if (name .eq. "i_th2_total_local") then
561  allocate(field_value%real_1d_array(current_state%local_grid%size(z_index)))
562  do k = 1, current_state%local_grid%size(z_index)
563  field_value%real_1d_array(k)=th2_tot(k)
564  enddo
565  end if
566  end subroutine field_value_retrieval_callback
567 end module profile_diagnostics_mod
real(kind=default_precision) qlcrit
integer, parameter, public component_scalar_field_type
Wrapper type for the value returned for a published field from a component.
real(kind=default_precision), dimension(:), allocatable uinit
real(kind=default_precision), dimension(:), allocatable ql_tot
type(standard_q_names_type), public standard_q_names
Definition: q_indices.F90:59
real(kind=default_precision), dimension(:), allocatable qv_tot
real(kind=default_precision), dimension(:), allocatable rh_tot
subroutine field_information_retrieval_callback(current_state, name, field_information)
Field information retrieval callback, this returns information for a specific components published fi...
real(kind=default_precision), dimension(:), allocatable rhon
real(kind=default_precision), dimension(:), allocatable uw_tot
integer, parameter, public default_precision
MPI communication type which we use for the prognostic and calculation data.
Definition: datadefn.F90:17
integer, parameter, public z_index
Grid index parameters.
Definition: grids.F90:14
subroutine timestep_callback(current_state)
real(kind=default_precision), dimension(:), allocatable thref
Contains common definitions for the data and datatypes used by MONC.
Definition: datadefn.F90:2
real(kind=default_precision), dimension(:), allocatable thinit
The ModelState which represents the current state of a run.
Definition: state.F90:39
real(kind=default_precision), dimension(:), allocatable tempfac
real(kind=default_precision), dimension(:), allocatable v_wind_tot
Saturation physics functionality which is used throughout the code.
Definition: saturation.F90:2
integer, parameter, public component_array_field_type
real(kind=default_precision), dimension(:), allocatable uprime_tot
real(kind=default_precision), dimension(:), allocatable vprime_tot
This manages the Q variables and specifically the mapping between names and the index that they are s...
Definition: q_indices.F90:2
Interfaces and types that MONC components must specify.
real(kind=default_precision), dimension(:), allocatable vinit
real(kind=default_precision), dimension(:), allocatable u_wind_tot
real(kind=default_precision), dimension(:), allocatable vw_tot
integer, parameter, public component_integer_data_type
real(kind=default_precision), dimension(:), allocatable theta_tot
real(kind=default_precision), dimension(:), allocatable rho
real(kind=default_precision), dimension(:), allocatable wq_tot
subroutine field_value_retrieval_callback(current_state, name, field_value)
Field value retrieval callback, this returns the value of a specific published field.
real(kind=default_precision) function, public qsaturation(temperature, pressure)
Function to return the saturation mixing ratio over water based on tetans formular QS=3...
Definition: saturation.F90:22
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.
real(kind=default_precision), dimension(:), allocatable prefn
Functionality to support the different types of grid and abstraction between global grids and local o...
Definition: grids.F90:5
type(component_descriptor_type) function, public profile_diagnostics_get_descriptor()
Provides the component descriptor for the core to register.
real(kind=default_precision), dimension(:), allocatable ww_tot
Manages the options database. Contains administration functions and deduce runtime options from the c...
real(kind=default_precision), dimension(:), allocatable wtheta_tot
real(kind=default_precision), dimension(:), allocatable th2_tot
real(kind=default_precision), dimension(:), allocatable w_wind_tot
The model state which represents the current state of a run.
Definition: state.F90:2
integer, parameter, public y_index
Definition: grids.F90:14
integer, parameter, public x_index
Definition: grids.F90:14
integer function, public get_q_index(name, assigning_component)
Add in a new entry into the register if the name does not already exist or return the index of the pr...
Definition: q_indices.F90:112
subroutine initialisation_callback(current_state)
integer, parameter, public component_double_data_type