MONC
stepfields.F90
Go to the documentation of this file.
1 
13  use science_constants_mod, only : rlargep
15 
16  implicit none
17 
18 #ifndef TEST_MODE
19  private
20 #endif
21 
23 
24 
25  real(kind=default_precision), allocatable :: resetq_min(:)
26  logical :: l_nonconservative_positive_q=.true.
27 
28  ! Local tendency diagnostic variables for this component
29  ! 3D tendency fields and logicals for their use
30  real(kind=default_precision), dimension(:,:,:), allocatable :: &
34  logical :: l_tend_3d_th,l_tend_3d_qv, &
37  ! Local mean tendency profile fields and logicals for their use
38  real(kind=default_precision), dimension(:), allocatable :: &
45  ! q indices
46  integer :: iqv=0, iql=0, iqr=0, iqi=0, iqs=0, iqg=0
48 
49 
51 
52 contains
53 
57  stepfields_get_descriptor%name="stepfields"
58  stepfields_get_descriptor%version=0.1
62 
65  allocate(stepfields_get_descriptor%published_fields(8+8))
66 
67  stepfields_get_descriptor%published_fields(1)="tend_th_total_3d_local"
68  stepfields_get_descriptor%published_fields(2)="tend_qv_total_3d_local"
69  stepfields_get_descriptor%published_fields(3)="tend_ql_total_3d_local"
70  stepfields_get_descriptor%published_fields(4)="tend_qi_total_3d_local"
71  stepfields_get_descriptor%published_fields(5)="tend_qr_total_3d_local"
72  stepfields_get_descriptor%published_fields(6)="tend_qs_total_3d_local"
73  stepfields_get_descriptor%published_fields(7)="tend_qg_total_3d_local"
74  stepfields_get_descriptor%published_fields(8)="tend_tabs_total_3d_local"
75 
76  stepfields_get_descriptor%published_fields(8+1)="tend_th_total_profile_total_local"
77  stepfields_get_descriptor%published_fields(8+2)="tend_qv_total_profile_total_local"
78  stepfields_get_descriptor%published_fields(8+3)="tend_ql_total_profile_total_local"
79  stepfields_get_descriptor%published_fields(8+4)="tend_qi_total_profile_total_local"
80  stepfields_get_descriptor%published_fields(8+5)="tend_qr_total_profile_total_local"
81  stepfields_get_descriptor%published_fields(8+6)="tend_qs_total_profile_total_local"
82  stepfields_get_descriptor%published_fields(8+7)="tend_qg_total_profile_total_local"
83  stepfields_get_descriptor%published_fields(8+8)="tend_tabs_total_profile_total_local"
84 
85  end function stepfields_get_descriptor
86 
89  subroutine initialisation_callback(current_state)
90  type(model_state_type), target, intent(inout) :: current_state
91  logical :: l_qdiag
92 
93 #ifdef U_ACTIVE
94  allocate(current_state%global_grid%configuration%vertical%savolubar(current_state%local_grid%size(z_index)))
95 #endif
96 #ifdef V_ACTIVE
97  allocate(current_state%global_grid%configuration%vertical%savolvbar(current_state%local_grid%size(z_index)))
98 #endif
99 
100  allocate(resetq_min(current_state%number_q_fields))
101  cfl_is_enabled=is_component_enabled(current_state%options_database, "cfltest")
102  if (cfl_is_enabled) call reset_local_minmax_values(current_state)
103 
104  ! Set tendency diagnostic logicals based on availability
105  ! Need to use 3d tendencies to compute the profiles, so they will be allocated
106  ! in the case where profiles are available
107  l_qdiag = (.not. current_state%passive_q .and. current_state%number_q_fields .gt. 0)
108 
109  l_tend_pr_tot_th = current_state%th%active
110  l_tend_pr_tot_qv = l_qdiag .and. current_state%number_q_fields .ge. 1
111  l_tend_pr_tot_ql = l_qdiag .and. current_state%number_q_fields .ge. 2
112  l_tend_pr_tot_qi = l_qdiag .and. current_state%number_q_fields .ge. 11
113  l_tend_pr_tot_qr = l_qdiag .and. current_state%number_q_fields .ge. 11
114  l_tend_pr_tot_qs = l_qdiag .and. current_state%number_q_fields .ge. 11
115  l_tend_pr_tot_qg = l_qdiag .and. current_state%number_q_fields .ge. 11
117 
118  l_tend_3d_th = current_state%th%active .or. l_tend_pr_tot_th
119  l_tend_3d_qv = (l_qdiag .and. current_state%number_q_fields .ge. 1) .or. l_tend_pr_tot_qv
120  l_tend_3d_ql = (l_qdiag .and. current_state%number_q_fields .ge. 2) .or. l_tend_pr_tot_ql
121  l_tend_3d_qi = (l_qdiag .and. current_state%number_q_fields .ge. 11) .or. l_tend_pr_tot_qi
122  l_tend_3d_qr = (l_qdiag .and. current_state%number_q_fields .ge. 11) .or. l_tend_pr_tot_qr
123  l_tend_3d_qs = (l_qdiag .and. current_state%number_q_fields .ge. 11) .or. l_tend_pr_tot_qs
124  l_tend_3d_qg = (l_qdiag .and. current_state%number_q_fields .ge. 11) .or. l_tend_pr_tot_qg
126 
127  ! Allocate 3d tendency fields upon availability
128  if (l_tend_3d_th) then
129  allocate( tend_3d_th(current_state%local_grid%size(z_index), &
130  current_state%local_grid%size(y_index), &
131  current_state%local_grid%size(x_index) ) )
132  endif
133  if (l_tend_3d_qv) then
134  iqv=get_q_index(standard_q_names%VAPOUR, 'stepfields')
135  allocate( tend_3d_qv(current_state%local_grid%size(z_index), &
136  current_state%local_grid%size(y_index), &
137  current_state%local_grid%size(x_index) ) )
138  endif
139  if (l_tend_3d_ql) then
140  iql=get_q_index(standard_q_names%CLOUD_LIQUID_MASS, 'stepfields')
141  allocate( tend_3d_ql(current_state%local_grid%size(z_index), &
142  current_state%local_grid%size(y_index), &
143  current_state%local_grid%size(x_index) ) )
144  endif
145  if (l_tend_3d_qi) then
146  iqi=get_q_index(standard_q_names%ICE_MASS, 'stepfields')
147  allocate( tend_3d_qi(current_state%local_grid%size(z_index), &
148  current_state%local_grid%size(y_index), &
149  current_state%local_grid%size(x_index) ) )
150  endif
151  if (l_tend_3d_qr) then
152  iqr=get_q_index(standard_q_names%RAIN_MASS, 'stepfields')
153  allocate( tend_3d_qr(current_state%local_grid%size(z_index), &
154  current_state%local_grid%size(y_index), &
155  current_state%local_grid%size(x_index) ) )
156  endif
157  if (l_tend_3d_qs) then
158  iqs=get_q_index(standard_q_names%SNOW_MASS, 'stepfields')
159  allocate( tend_3d_qs(current_state%local_grid%size(z_index), &
160  current_state%local_grid%size(y_index), &
161  current_state%local_grid%size(x_index) ) )
162  endif
163  if (l_tend_3d_qg) then
164  iqg=get_q_index(standard_q_names%GRAUPEL_MASS, 'stepfields')
165  allocate( tend_3d_qg(current_state%local_grid%size(z_index), &
166  current_state%local_grid%size(y_index), &
167  current_state%local_grid%size(x_index) ) )
168  endif
169  if (l_tend_3d_tabs) then
170  allocate( tend_3d_tabs(current_state%local_grid%size(z_index), &
171  current_state%local_grid%size(y_index), &
172  current_state%local_grid%size(x_index) ) )
173  endif
174 
175  ! Allocate profile tendency fields upon availability
176  if (l_tend_pr_tot_th) then
177  allocate( tend_pr_tot_th(current_state%local_grid%size(z_index)) )
178  endif
179  if (l_tend_pr_tot_qv) then
180  allocate( tend_pr_tot_qv(current_state%local_grid%size(z_index)) )
181  endif
182  if (l_tend_pr_tot_ql) then
183  allocate( tend_pr_tot_ql(current_state%local_grid%size(z_index)) )
184  endif
185  if (l_tend_pr_tot_qi) then
186  allocate( tend_pr_tot_qi(current_state%local_grid%size(z_index)) )
187  endif
188  if (l_tend_pr_tot_qr) then
189  allocate( tend_pr_tot_qr(current_state%local_grid%size(z_index)) )
190  endif
191  if (l_tend_pr_tot_qs) then
192  allocate( tend_pr_tot_qs(current_state%local_grid%size(z_index)) )
193  endif
194  if (l_tend_pr_tot_qg) then
195  allocate( tend_pr_tot_qg(current_state%local_grid%size(z_index)) )
196  endif
197  if (l_tend_pr_tot_tabs) then
198  allocate( tend_pr_tot_tabs(current_state%local_grid%size(z_index)) )
199  endif
200 
201  ! Save the sampling_frequency to force diagnostic calculation on select time steps
202  diagnostic_generation_frequency=options_get_integer(current_state%options_database, "sampling_frequency")
203 
204  end subroutine initialisation_callback
205 
206 
209  subroutine finalisation_callback(current_state)
210  type(model_state_type), target, intent(inout) :: current_state
211 
212  deallocate(resetq_min)
213 
214  if (allocated(tend_3d_th)) deallocate(tend_3d_th)
215  if (allocated(tend_3d_qv)) deallocate(tend_3d_qv)
216  if (allocated(tend_3d_ql)) deallocate(tend_3d_ql)
217  if (allocated(tend_3d_qi)) deallocate(tend_3d_qi)
218  if (allocated(tend_3d_qr)) deallocate(tend_3d_qr)
219  if (allocated(tend_3d_qs)) deallocate(tend_3d_qs)
220  if (allocated(tend_3d_qg)) deallocate(tend_3d_qg)
221  if (allocated(tend_3d_tabs)) deallocate(tend_3d_tabs)
222 
223  if (allocated(tend_pr_tot_th)) deallocate(tend_pr_tot_th)
224  if (allocated(tend_pr_tot_qv)) deallocate(tend_pr_tot_qv)
225  if (allocated(tend_pr_tot_ql)) deallocate(tend_pr_tot_ql)
226  if (allocated(tend_pr_tot_qi)) deallocate(tend_pr_tot_qi)
227  if (allocated(tend_pr_tot_qr)) deallocate(tend_pr_tot_qr)
228  if (allocated(tend_pr_tot_qs)) deallocate(tend_pr_tot_qs)
229  if (allocated(tend_pr_tot_qg)) deallocate(tend_pr_tot_qg)
230  if (allocated(tend_pr_tot_tabs)) deallocate(tend_pr_tot_tabs)
231 
232  end subroutine finalisation_callback
233 
236  subroutine timestep_callback(current_state)
237  type(model_state_type), target, intent(inout) :: current_state
238 
239  integer :: iq
240  integer :: current_x_index, current_y_index, target_x_index, target_y_index
241 
242  current_x_index=current_state%column_local_x
243  current_y_index=current_state%column_local_y
244  target_y_index=current_y_index-current_state%local_grid%halo_size(y_index)
245  target_x_index=current_x_index-current_state%local_grid%halo_size(x_index)
246 
247 
248  if (cfl_is_enabled .and. current_state%first_timestep_column) then
249  if (mod(current_state%timestep, current_state%cfl_frequency) == 1 .or. &
250  current_state%timestep-current_state%start_timestep .le. current_state%cfl_frequency) then
251  determine_flow_minmax=.true.
252  call reset_local_minmax_values(current_state)
253  else
254  determine_flow_minmax=.false.
255  end if
256  end if
257 
258  if (.not. current_state%halo_column) then
260  call determine_local_flow_minmax(current_state, current_state%column_local_y, current_state%column_local_x)
261  call step_all_fields(current_state)
262  end if
263 
264  ! Remove negative rounding errors
266  do iq=1,current_state%number_q_fields
267  if (current_state%first_timestep_column)then
268  resetq_min(iq)=minval(current_state%zq(iq)%data(:,current_state%column_local_y, current_state%column_local_x))
269  else
270  resetq_min(iq)=min(resetq_min(iq),&
271  minval(current_state%zq(iq)%data(:,current_state%column_local_y, current_state%column_local_x)))
272  end if
273  call remove_negative_rounding_errors_for_single_field(current_state%column_local_x, current_state%column_local_y, &
274  current_state%column_local_x-2, current_state%column_local_y-1, current_state%zq(iq), current_state%local_grid)
275  end do
276  end if
277 
278 
279  ! Zero profile tendency totals on first instance in the sum
280  if (current_state%first_timestep_column) then
281  if (l_tend_pr_tot_th) then
282  tend_pr_tot_th(:)=0.0_default_precision
283  endif
284  if (l_tend_pr_tot_qv) then
285  tend_pr_tot_qv(:)=0.0_default_precision
286  endif
287  if (l_tend_pr_tot_ql) then
288  tend_pr_tot_ql(:)=0.0_default_precision
289  endif
290  if (l_tend_pr_tot_qi) then
291  tend_pr_tot_qi(:)=0.0_default_precision
292  endif
293  if (l_tend_pr_tot_qr) then
294  tend_pr_tot_qr(:)=0.0_default_precision
295  endif
296  if (l_tend_pr_tot_qs) then
297  tend_pr_tot_qs(:)=0.0_default_precision
298  endif
299  if (l_tend_pr_tot_qg) then
300  tend_pr_tot_qg(:)=0.0_default_precision
301  endif
302  if (l_tend_pr_tot_tabs) then
303  tend_pr_tot_tabs(:)=0.0_default_precision
304  endif
305  endif ! zero totals
306 
307  if (mod(current_state%timestep, diagnostic_generation_frequency) == 0 .and. .not. current_state%halo_column) then
308  call compute_component_tendencies(current_state, current_x_index, current_y_index, target_x_index, target_y_index)
309  end if
310 
311  end subroutine timestep_callback
312 
315  subroutine step_all_fields(current_state)
316  type(model_state_type), target, intent(inout) :: current_state
317 
318  integer :: x_prev, y_prev, i, k
319  real(kind=default_precision) :: c1, c2
320 
321  x_prev = current_state%column_local_x-2
322  y_prev = current_state%column_local_y-1
323 
324  c1 = 1.0_default_precision - 2.0_default_precision*current_state%tsmth
325  c2 = current_state%tsmth
326 
327 #ifdef U_ACTIVE
328  current_state%global_grid%configuration%vertical%savolubar=current_state%global_grid%configuration%vertical%olzubar
329  call step_single_field(current_state%column_local_x, current_state%column_local_y, &
330  x_prev, y_prev, current_state%u, current_state%zu, current_state%su, current_state%local_grid, .true., &
331  current_state%field_stepping, current_state%dtm, current_state%ugal, c1, c2, .false., current_state%savu)
332 #endif
333 #ifdef V_ACTIVE
334  current_state%global_grid%configuration%vertical%savolvbar=current_state%global_grid%configuration%vertical%olzvbar
335  call step_single_field(current_state%column_local_x, current_state%column_local_y, &
336  x_prev, y_prev, current_state%v, current_state%zv, current_state%sv, current_state%local_grid, .true., &
337  current_state%field_stepping, current_state%dtm, current_state%vgal, c1, c2, .false., current_state%savv)
338 #endif
339 #ifdef W_ACTIVE
340  call step_single_field(current_state%column_local_x, current_state%column_local_y, &
341  x_prev, y_prev, current_state%w, current_state%zw, current_state%sw, current_state%local_grid, .false., &
342  current_state%field_stepping, current_state%dtm, real(0., kind=default_precision), c1, c2, .false., current_state%savw)
343 #endif
344  if (current_state%th%active) then
345  call step_single_field(current_state%column_local_x, current_state%column_local_y, &
346  x_prev, y_prev, current_state%th, current_state%zth, current_state%sth, current_state%local_grid, .false., &
347  current_state%field_stepping, current_state%dtm, real(0., kind=default_precision), c1, c2, &
348  current_state%field_stepping == centred_stepping)
349  endif
350  do i=1,current_state%number_q_fields
351  if (current_state%q(i)%active) then
352  call step_single_field(current_state%column_local_x, current_state%column_local_y, x_prev, y_prev, &
353  current_state%q(i), current_state%zq(i), current_state%sq(i), current_state%local_grid, .false., &
354  current_state%field_stepping, current_state%dtm, real(0., kind=default_precision), c1, c2, &
355  current_state%field_stepping == centred_stepping)
356  end if
357  end do
358  end subroutine step_all_fields
359 
365  subroutine determine_local_flow_minmax(current_state, local_y, local_x)
366  type(model_state_type), target, intent(inout) :: current_state
367  integer, intent(in) :: local_x, local_y
368 
369  integer :: k
370 
371  do k=2, current_state%local_grid%local_domain_end_index(z_index)
372 #ifdef U_ACTIVE
373  current_state%local_zumax = max(current_state%local_zumax, current_state%zu%data(k,local_y,local_x))
374  current_state%local_zumin = min(current_state%local_zumin, current_state%zu%data(k,local_y,local_x))
375 #endif
376 #ifdef V_ACTIVE
377  current_state%local_zvmax = max(current_state%local_zvmax, current_state%zv%data(k,local_y,local_x))
378  current_state%local_zvmin = min(current_state%local_zvmin, current_state%zv%data(k,local_y,local_x))
379 #endif
380 #ifdef W_ACTIVE
381  if (k .lt. current_state%local_grid%local_domain_end_index(z_index)) then
382  current_state%abswmax(k) = max(current_state%abswmax(k), abs(current_state%zw%data(k,local_y,local_x)))
383  end if
384 #endif
385  end do
386  end subroutine determine_local_flow_minmax
387 
390  subroutine reset_local_minmax_values(current_state)
391  type(model_state_type), intent(inout), target :: current_state
392 
393  ! Reset the local values for the next timestep
394  current_state%local_zumin=rlargep
395  current_state%local_zumax=-rlargep
396  current_state%local_zvmin=rlargep
397  current_state%local_zvmax=-rlargep
398  current_state%abswmax=-rlargep
399  end subroutine reset_local_minmax_values
400 
409  subroutine remove_negative_rounding_errors_for_single_field(x_local_index, y_local_index, x_prev, y_prev, field, local_grid)
410  integer, intent(in) :: x_local_index, y_local_index, x_prev, y_prev
411  type(local_grid_type), intent(inout) :: local_grid
412  type(prognostic_field_type), intent(inout) :: field
413 
414  if (x_prev .ge. local_grid%local_domain_start_index(x_index)) then
415  call remove_negative_rounding_errors_in_slice(y_local_index, x_prev, y_prev, field, local_grid)
416  end if
417  if (x_local_index == local_grid%local_domain_end_index(x_index)) then
418  if (x_local_index .gt. 1) then
419  call remove_negative_rounding_errors_in_slice(y_local_index, x_local_index-1, y_prev, field, local_grid)
420  end if
421  call remove_negative_rounding_errors_in_slice(y_local_index, x_local_index, y_prev, field, local_grid)
422  end if
424 
432  subroutine remove_negative_rounding_errors_in_slice(y_local_index, x_prev, y_prev, field, local_grid)
433  integer, intent(in) :: y_local_index, x_prev, y_prev
434  type(local_grid_type), intent(inout) :: local_grid
435  type(prognostic_field_type), intent(inout) :: field
436 
437  if (y_prev .ge. local_grid%local_domain_start_index(y_index)) then
438  where (field%data(:, y_prev, x_prev) < 0.0_default_precision)
439  field%data(:, y_prev, x_prev)=0.0_default_precision
440  end where
441  end if
442  if (y_local_index == local_grid%local_domain_end_index(y_index)) then
443  where (field%data(:, y_local_index, x_prev) < 0.0_default_precision)
444  field%data(:, y_local_index, x_prev)=0.0_default_precision
445  end where
446  end if
448 
466  subroutine step_single_field(x_local_index, y_local_index, x_prev, y_prev, field, zfield, sfield, local_grid,&
467  flow_field, direction, dtm, gal, c1, c2, do_timesmoothing, sav)
468  integer, intent(in) :: x_local_index, y_local_index, x_prev, y_prev, direction
469  real(kind=default_precision), intent(in) :: dtm, gal
470  logical, intent(in) :: flow_field, do_timesmoothing
471  type(local_grid_type), intent(inout) :: local_grid
472  type(prognostic_field_type), intent(inout) :: field, zfield, sfield
473  real(kind=default_precision), intent(in) :: c1, c2
474  type(prognostic_field_type), optional, intent(inout) :: sav
475 
476  if (x_prev .ge. local_grid%local_domain_start_index(x_index)) then
477  if (present(sav)) then
478  call step_column_in_slice(y_local_index, x_prev, y_prev, field, zfield, sfield, local_grid, &
479  flow_field, direction, dtm, gal, c1, c2, do_timesmoothing, sav)
480  else
481  call step_column_in_slice(y_local_index, x_prev, y_prev, field, zfield, sfield, local_grid, &
482  flow_field, direction, dtm, gal, c1, c2, do_timesmoothing)
483  end if
484  end if
485 
486  if (x_local_index == local_grid%local_domain_end_index(x_index)) then
487  ! If this is the last slice then process x-1 (if applicable) and x too
488  if (x_local_index .gt. 1) then
489  if (present(sav)) then
490  call step_column_in_slice(y_local_index, x_local_index-1, y_prev, field, zfield, sfield, local_grid, &
491  flow_field, direction, dtm, gal, c1, c2, do_timesmoothing, sav)
492  else
493  call step_column_in_slice(y_local_index, x_local_index-1, y_prev, field, zfield, sfield, local_grid, &
494  flow_field, direction, dtm, gal, c1, c2, do_timesmoothing)
495  end if
496  end if
497  if (present(sav)) then
498  call step_column_in_slice(y_local_index, x_local_index, y_prev, field, zfield, sfield, local_grid, &
499  flow_field, direction, dtm, gal, c1, c2, do_timesmoothing, sav)
500  else
501  call step_column_in_slice(y_local_index, x_local_index, y_prev, field, zfield, sfield, local_grid, &
502  flow_field, direction, dtm, gal, c1, c2, do_timesmoothing)
503  end if
504  end if
505  end subroutine step_single_field
506 
515  subroutine perform_timesmooth_for_field(field, zfield, local_grid, x_index, y_index, c1, c2)
516  type(prognostic_field_type), intent(inout) :: field, zfield
517  type(local_grid_type), intent(inout) :: local_grid
518  integer, intent(in) :: x_index, y_index
519  real(kind=default_precision), intent(in) :: c1, c2
520 
521  integer :: k
522 
523  do k=1,local_grid%size(z_index)
524  field%data(k, y_index, x_index)=c1*field%data(k, y_index, x_index)+c2*zfield%data(k, y_index, x_index)
525  end do
526  end subroutine perform_timesmooth_for_field
527 
546  subroutine step_column_in_slice(y_local_index, x_prev, y_prev, field, zfield, sfield, local_grid,&
547  flow_field, direction, dtm, gal, c1, c2, do_timesmoothing, sav)
548  integer, intent(in) :: y_local_index, x_prev, y_prev, direction
549  real(kind=default_precision), intent(in) :: dtm, gal
550  logical, intent(in) :: flow_field, do_timesmoothing
551  type(local_grid_type), intent(inout) :: local_grid
552  type(prognostic_field_type), intent(inout) :: field, zfield, sfield
553  real(kind=default_precision), intent(in) :: c1, c2
554  type(prognostic_field_type), optional, intent(inout) :: sav
555 
556  if (y_prev .ge. local_grid%local_domain_start_index(y_index)) then
557  if (do_timesmoothing) then
558  call perform_timesmooth_for_field(field, zfield, local_grid, x_prev, y_prev, c1, c2)
559  end if
560  if (present(sav)) then
561  call step_field(x_prev, y_prev, field, zfield, sfield, local_grid, flow_field, direction, dtm, gal, sav)
562  else
563  call step_field(x_prev, y_prev, field, zfield, sfield, local_grid, flow_field, direction, dtm, gal)
564  end if
565  end if
566 
567  if (y_local_index == local_grid%local_domain_end_index(y_index)) then
568  if (do_timesmoothing) then
569  call perform_timesmooth_for_field(field, zfield, local_grid, x_prev, y_local_index, c1, c2)
570  end if
571  if (present(sav)) then
572  call step_field(x_prev, y_local_index, field, zfield, sfield, local_grid, flow_field, direction, dtm, gal, sav)
573  else
574  call step_field(x_prev, y_local_index, field, zfield, sfield, local_grid, flow_field, direction, dtm, gal)
575  end if
576  end if
577  end subroutine step_column_in_slice
578 
590  subroutine step_field(x_local_index, y_local_index, field, zfield, sfield, local_grid, flow_field, direction, dtm, gal, sav)
591  integer, intent(in) :: x_local_index, y_local_index, direction
592  real(kind=default_precision), intent(in) :: dtm, gal
593  logical, intent(in) :: flow_field
594  type(local_grid_type), intent(inout) :: local_grid
595  type(prognostic_field_type), intent(inout) :: field, zfield, sfield
596  type(prognostic_field_type), optional, intent(inout) :: sav
597 
598  integer :: k
599  real(kind=default_precision) :: actual_gal, dtm_x2
600 
601  dtm_x2 = 2.0_default_precision * dtm
602 
603  actual_gal = merge(gal, real(0.0_default_precision, kind=default_precision), flow_field)
604 
605  sfield%data(1,y_local_index, x_local_index)=0.0_default_precision
606 
607  do k=1,local_grid%size(z_index)
608  ! Save the Z field which is used in the Robert filter
609  if (present(sav) .and. direction .eq. centred_stepping) &
610  sav%data(k,y_local_index, x_local_index) = zfield%data(k, y_local_index, x_local_index) + actual_gal
611  if (flow_field) field%data(k, y_local_index, x_local_index) = actual_gal + field%data(k, y_local_index, x_local_index)
612  if (direction == forward_stepping) then
613  zfield%data(k, y_local_index, x_local_index) = field%data(k, y_local_index, x_local_index) + dtm * &
614  sfield%data(k, y_local_index, x_local_index)
615  else
616  zfield%data(k, y_local_index, x_local_index) = actual_gal+zfield%data(k, y_local_index, x_local_index)+dtm_x2*&
617  sfield%data(k, y_local_index, x_local_index)
618  end if
619  end do
620  end subroutine step_field
621 
628  subroutine compute_component_tendencies(current_state, cxn, cyn, txn, tyn)
629  type(model_state_type), target, intent(inout) :: current_state
630  integer, intent(in) :: cxn, cyn, txn, tyn
631 
632  ! Calculate change in tendency due to component
633  if (l_tend_3d_th) then
634  tend_3d_th(:,tyn,txn)=current_state%sth%data(:,cyn,cxn)
635  endif
636  if (l_tend_3d_qv) then
637  tend_3d_qv(:,tyn,txn)=current_state%sq(iqv)%data(:,cyn,cxn)
638  endif
639  if (l_tend_3d_ql) then
640  tend_3d_ql(:,tyn,txn)=current_state%sq(iql)%data(:,cyn,cxn)
641  endif
642  if (l_tend_3d_qi) then
643  tend_3d_qi(:,tyn,txn)=current_state%sq(iqi)%data(:,cyn,cxn)
644  endif
645  if (l_tend_3d_qr) then
646  tend_3d_qr(:,tyn,txn)=current_state%sq(iqr)%data(:,cyn,cxn)
647  endif
648  if (l_tend_3d_qs) then
649  tend_3d_qs(:,tyn,txn)=current_state%sq(iqs)%data(:,cyn,cxn)
650  endif
651  if (l_tend_3d_qg) then
652  tend_3d_qg(:,tyn,txn)=current_state%sq(iqg)%data(:,cyn,cxn)
653  endif
654  if (l_tend_3d_tabs) then
655  tend_3d_tabs(:,tyn,txn)= &
656  current_state%sth%data(:,cyn,cxn) * current_state%global_grid%configuration%vertical%rprefrcp(:)
657  endif
658 
659  ! Add local tendency fields to the profile total
660  if (l_tend_pr_tot_th) then
661  tend_pr_tot_th(:)=tend_pr_tot_th(:) + tend_3d_th(:,tyn,txn)
662  endif
663  if (l_tend_pr_tot_qv) then
664  tend_pr_tot_qv(:)=tend_pr_tot_qv(:) + tend_3d_qv(:,tyn,txn)
665  endif
666  if (l_tend_pr_tot_ql) then
667  tend_pr_tot_ql(:)=tend_pr_tot_ql(:) + tend_3d_ql(:,tyn,txn)
668  endif
669  if (l_tend_pr_tot_qi) then
670  tend_pr_tot_qi(:)=tend_pr_tot_qi(:) + tend_3d_qi(:,tyn,txn)
671  endif
672  if (l_tend_pr_tot_qr) then
673  tend_pr_tot_qr(:)=tend_pr_tot_qr(:) + tend_3d_qr(:,tyn,txn)
674  endif
675  if (l_tend_pr_tot_qs) then
676  tend_pr_tot_qs(:)=tend_pr_tot_qs(:) + tend_3d_qs(:,tyn,txn)
677  endif
678  if (l_tend_pr_tot_qg) then
679  tend_pr_tot_qg(:)=tend_pr_tot_qg(:) + tend_3d_qg(:,tyn,txn)
680  endif
681  if (l_tend_pr_tot_tabs) then
683  endif
684 
685  end subroutine compute_component_tendencies
686 
687 
693  subroutine field_information_retrieval_callback(current_state, name, field_information)
694  type(model_state_type), target, intent(inout) :: current_state
695  character(len=*), intent(in) :: name
696  type(component_field_information_type), intent(out) :: field_information
697  integer :: strcomp
698 
699  ! Field information for 3d
700  strcomp=index(name, "_total_3d_local")
701  if (strcomp .ne. 0) then
702  field_information%field_type=component_array_field_type
703  field_information%number_dimensions=3
704  field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
705  field_information%dimension_sizes(2)=current_state%local_grid%size(y_index)
706  field_information%dimension_sizes(3)=current_state%local_grid%size(x_index)
707  field_information%data_type=component_double_data_type
708 
709  if (name .eq. "tend_th_total_3d_local") then
710  field_information%enabled=l_tend_3d_th
711  else if (name .eq. "tend_qv_total_3d_local") then
712  field_information%enabled=l_tend_3d_qv
713  else if (name .eq. "tend_ql_total_3d_local") then
714  field_information%enabled=l_tend_3d_ql
715  else if (name .eq. "tend_qi_total_3d_local") then
716  field_information%enabled=l_tend_3d_qi
717  else if (name .eq. "tend_qr_total_3d_local") then
718  field_information%enabled=l_tend_3d_qr
719  else if (name .eq. "tend_qs_total_3d_local") then
720  field_information%enabled=l_tend_3d_qs
721  else if (name .eq. "tend_qg_total_3d_local") then
722  field_information%enabled=l_tend_3d_qg
723  else if (name .eq. "tend_tabs_total_3d_local") then
724  field_information%enabled=l_tend_3d_tabs
725  else
726  field_information%enabled=.true.
727  end if
728 
729  end if !end 3d check
730 
731  ! Field information for profiles
732  strcomp=index(name, "_total_profile_total_local")
733  if (strcomp .ne. 0) then
734  field_information%field_type=component_array_field_type
735  field_information%number_dimensions=1
736  field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
737  field_information%data_type=component_double_data_type
738 
739  if (name .eq. "tend_th_total_profile_total_local") then
740  field_information%enabled=l_tend_pr_tot_th
741  else if (name .eq. "tend_qv_total_profile_total_local") then
742  field_information%enabled=l_tend_pr_tot_qv
743  else if (name .eq. "tend_ql_total_profile_total_local") then
744  field_information%enabled=l_tend_pr_tot_ql
745  else if (name .eq. "tend_qi_total_profile_total_local") then
746  field_information%enabled=l_tend_pr_tot_qi
747  else if (name .eq. "tend_qr_total_profile_total_local") then
748  field_information%enabled=l_tend_pr_tot_qr
749  else if (name .eq. "tend_qs_total_profile_total_local") then
750  field_information%enabled=l_tend_pr_tot_qs
751  else if (name .eq. "tend_qg_total_profile_total_local") then
752  field_information%enabled=l_tend_pr_tot_qg
753  else if (name .eq. "tend_tabs_total_profile_total_local") then
754  field_information%enabled=l_tend_pr_tot_tabs
755  else
756  field_information%enabled=.true.
757  end if
758 
759  end if !end profile check
760 
762 
763 
768  subroutine field_value_retrieval_callback(current_state, name, field_value)
769  type(model_state_type), target, intent(inout) :: current_state
770  character(len=*), intent(in) :: name
771  type(component_field_value_type), intent(out) :: field_value
772 
773  ! 3d Tendency Fields
774  if (name .eq. "tend_th_total_3d_local" .and. allocated(tend_3d_th)) then
775  call set_published_field_value(field_value, real_3d_field=tend_3d_th)
776  else if (name .eq. "tend_qv_total_3d_local" .and. allocated(tend_3d_qv)) then
777  call set_published_field_value(field_value, real_3d_field=tend_3d_qv)
778  else if (name .eq. "tend_ql_total_3d_local" .and. allocated(tend_3d_ql)) then
779  call set_published_field_value(field_value, real_3d_field=tend_3d_ql)
780  else if (name .eq. "tend_qi_total_3d_local" .and. allocated(tend_3d_qi)) then
781  call set_published_field_value(field_value, real_3d_field=tend_3d_qi)
782  else if (name .eq. "tend_qr_total_3d_local" .and. allocated(tend_3d_qr)) then
783  call set_published_field_value(field_value, real_3d_field=tend_3d_qr)
784  else if (name .eq. "tend_qs_total_3d_local" .and. allocated(tend_3d_qs)) then
785  call set_published_field_value(field_value, real_3d_field=tend_3d_qs)
786  else if (name .eq. "tend_qg_total_3d_local" .and. allocated(tend_3d_qg)) then
787  call set_published_field_value(field_value, real_3d_field=tend_3d_qg)
788  else if (name .eq. "tend_tabs_total_3d_local" .and. allocated(tend_3d_tabs)) then
789  call set_published_field_value(field_value, real_3d_field=tend_3d_tabs)
790 
791  ! Profile Tendency Fields
792  else if (name .eq. "tend_th_total_profile_total_local" .and. allocated(tend_pr_tot_th)) then
793  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_th)
794  else if (name .eq. "tend_qv_total_profile_total_local" .and. allocated(tend_pr_tot_qv)) then
795  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_qv)
796  else if (name .eq. "tend_ql_total_profile_total_local" .and. allocated(tend_pr_tot_ql)) then
797  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_ql)
798  else if (name .eq. "tend_qi_total_profile_total_local" .and. allocated(tend_pr_tot_qi)) then
799  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_qi)
800  else if (name .eq. "tend_qr_total_profile_total_local" .and. allocated(tend_pr_tot_qr)) then
801  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_qr)
802  else if (name .eq. "tend_qs_total_profile_total_local" .and. allocated(tend_pr_tot_qs)) then
803  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_qs)
804  else if (name .eq. "tend_qg_total_profile_total_local" .and. allocated(tend_pr_tot_qg)) then
805  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_qg)
806  else if (name .eq. "tend_tabs_total_profile_total_local" .and. allocated(tend_pr_tot_tabs)) then
807  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_tabs)
808  end if
809 
810  end subroutine field_value_retrieval_callback
811 
812 
817  subroutine set_published_field_value(field_value, real_1d_field, real_2d_field, real_3d_field)
818  type(component_field_value_type), intent(inout) :: field_value
819  real(kind=default_precision), dimension(:), optional :: real_1d_field
820  real(kind=default_precision), dimension(:,:), optional :: real_2d_field
821  real(kind=default_precision), dimension(:,:,:), optional :: real_3d_field
822 
823  if (present(real_1d_field)) then
824  allocate(field_value%real_1d_array(size(real_1d_field)), source=real_1d_field)
825  else if (present(real_2d_field)) then
826  allocate(field_value%real_2d_array(size(real_2d_field, 1), size(real_2d_field, 2)), source=real_2d_field)
827  else if (present(real_3d_field)) then
828  allocate(field_value%real_3d_array(size(real_3d_field, 1), size(real_3d_field, 2), size(real_3d_field, 3)), &
829  source=real_3d_field)
830  end if
831  end subroutine set_published_field_value
832 
833 end module stepfields_mod
registry_mod::is_component_enabled
logical function, public is_component_enabled(options_database, component_name)
Determines whether or not a specific component is registered and enabled.
Definition: registry.F90:334
science_constants_mod::rlargep
real(kind=default_precision), public rlargep
Definition: scienceconstants.F90:13
stepfields_mod::determine_flow_minmax
logical determine_flow_minmax
Definition: stepfields.F90:22
stepfields_mod::l_tend_pr_tot_qg
logical l_tend_pr_tot_qg
Definition: stepfields.F90:42
prognostics_mod
Contains prognostic field definitions and functions.
Definition: prognostics.F90:2
stepfields_mod::remove_negative_rounding_errors_for_single_field
subroutine remove_negative_rounding_errors_for_single_field(x_local_index, y_local_index, x_prev, y_prev, field, local_grid)
Removes the negative rounding errors from a specific single field. This works two columns behind and ...
Definition: stepfields.F90:410
stepfields_mod::l_tend_pr_tot_qs
logical l_tend_pr_tot_qs
Definition: stepfields.F90:42
stepfields_mod::iqv
integer iqv
Definition: stepfields.F90:46
stepfields_mod::l_tend_3d_tabs
logical l_tend_3d_tabs
Definition: stepfields.F90:34
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
stepfields_mod::l_tend_3d_ql
logical l_tend_3d_ql
Definition: stepfields.F90:34
grids_mod::y_index
integer, parameter, public y_index
Definition: grids.F90:14
stepfields_mod::stepfields_get_descriptor
type(component_descriptor_type) function, public stepfields_get_descriptor()
Provides the descriptor back to the caller and is used in component registration.
Definition: stepfields.F90:57
stepfields_mod::step_field
subroutine step_field(x_local_index, y_local_index, field, zfield, sfield, local_grid, flow_field, direction, dtm, gal, sav)
Will do the actual field stepping.
Definition: stepfields.F90:591
stepfields_mod::l_tend_3d_qr
logical l_tend_3d_qr
Definition: stepfields.F90:34
stepfields_mod::tend_3d_qg
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_qg
Definition: stepfields.F90:30
stepfields_mod::compute_component_tendencies
subroutine compute_component_tendencies(current_state, cxn, cyn, txn, tyn)
Computation of component tendencies.
Definition: stepfields.F90:629
stepfields_mod::tend_3d_qs
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_qs
Definition: stepfields.F90:30
stepfields_mod::tend_3d_qi
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_qi
Definition: stepfields.F90:30
stepfields_mod::tend_pr_tot_th
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_th
Definition: stepfields.F90:38
stepfields_mod::l_tend_3d_qv
logical l_tend_3d_qv
Definition: stepfields.F90:34
stepfields_mod::tend_3d_qv
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_qv
Definition: stepfields.F90:30
stepfields_mod::set_published_field_value
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.
Definition: stepfields.F90:818
stepfields_mod::l_tend_3d_th
logical l_tend_3d_th
Definition: stepfields.F90:34
stepfields_mod::cfl_is_enabled
logical cfl_is_enabled
Definition: stepfields.F90:22
stepfields_mod::step_all_fields
subroutine step_all_fields(current_state)
Steps all fields.
Definition: stepfields.F90:316
stepfields_mod::determine_local_flow_minmax
subroutine determine_local_flow_minmax(current_state, local_y, local_x)
Determines the minimum and maximum values for the local flow field. These are before the stepping,...
Definition: stepfields.F90:366
stepfields_mod::l_tend_3d_qs
logical l_tend_3d_qs
Definition: stepfields.F90:34
stepfields_mod::iqs
integer iqs
Definition: stepfields.F90:46
stepfields_mod::l_tend_3d_qi
logical l_tend_3d_qi
Definition: stepfields.F90:34
monc_component_mod
Interfaces and types that MONC components must specify.
Definition: monc_component.F90:6
stepfields_mod::field_value_retrieval_callback
subroutine field_value_retrieval_callback(current_state, name, field_value)
Field value retrieval callback, this returns the value of a specific published field.
Definition: stepfields.F90:769
stepfields_mod::tend_3d_ql
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_ql
Definition: stepfields.F90:30
stepfields_mod::step_column_in_slice
subroutine step_column_in_slice(y_local_index, x_prev, y_prev, field, zfield, sfield, local_grid, flow_field, direction, dtm, gal, c1, c2, do_timesmoothing, sav)
Will step a column in a specific slice. If y_prev is large enough then will step the y-1 column and i...
Definition: stepfields.F90:548
stepfields_mod::l_tend_pr_tot_tabs
logical l_tend_pr_tot_tabs
Definition: stepfields.F90:42
stepfields_mod::reset_local_minmax_values
subroutine reset_local_minmax_values(current_state)
Resets the local min and max values for the flow fields.
Definition: stepfields.F90:391
stepfields_mod::l_tend_pr_tot_qr
logical l_tend_pr_tot_qr
Definition: stepfields.F90:42
stepfields_mod::iqr
integer iqr
Definition: stepfields.F90:46
stepfields_mod::tend_pr_tot_qi
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_qi
Definition: stepfields.F90:38
stepfields_mod::l_tend_pr_tot_ql
logical l_tend_pr_tot_ql
Definition: stepfields.F90:42
stepfields_mod::iql
integer iql
Definition: stepfields.F90:46
stepfields_mod::step_single_field
subroutine step_single_field(x_local_index, y_local_index, x_prev, y_prev, field, zfield, sfield, local_grid, flow_field, direction, dtm, gal, c1, c2, do_timesmoothing, sav)
Steps a single specific field. This will step on the yth column of the x-2 slice and x-1 and x if thi...
Definition: stepfields.F90:468
stepfields_mod::finalisation_callback
subroutine finalisation_callback(current_state)
Finalisation callback.
Definition: stepfields.F90:210
science_constants_mod
Scientific constant values used throughout simulations. Each has a default value and this can be over...
Definition: scienceconstants.F90:3
grids_mod::z_index
integer, parameter, public z_index
Grid index parameters.
Definition: grids.F90:14
stepfields_mod::tend_pr_tot_ql
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_ql
Definition: stepfields.F90:38
stepfields_mod::l_tend_pr_tot_qv
logical l_tend_pr_tot_qv
Definition: stepfields.F90:42
q_indices_mod::standard_q_names
type(standard_q_names_type), public standard_q_names
Definition: q_indices.F90:59
stepfields_mod::tend_pr_tot_qr
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_qr
Definition: stepfields.F90:38
stepfields_mod::l_tend_pr_tot_qi
logical l_tend_pr_tot_qi
Definition: stepfields.F90:42
monc_component_mod::component_field_value_type
Wrapper type for the value returned for a published field from a component.
Definition: monc_component.F90:22
state_mod::model_state_type
The ModelState which represents the current state of a run.
Definition: state.F90:39
state_mod::forward_stepping
integer, parameter, public forward_stepping
Definition: state.F90:15
monc_component_mod::component_field_information_type
Definition: monc_component.F90:31
stepfields_mod::tend_pr_tot_qg
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_qg
Definition: stepfields.F90:38
stepfields_mod::tend_pr_tot_qv
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_qv
Definition: stepfields.F90:38
state_mod::centred_stepping
integer, parameter, public centred_stepping
Stepping parameter values which determine centred or forward stepping.
Definition: state.F90:15
grids_mod::local_grid_type
Defined the local grid, i.e. the grid held on this process after decomposition.
Definition: grids.F90:118
stepfields_mod::resetq_min
real(kind=default_precision), dimension(:), allocatable resetq_min
Definition: stepfields.F90:25
stepfields_mod::tend_3d_qr
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_qr
Definition: stepfields.F90:30
stepfields_mod::remove_negative_rounding_errors_in_slice
subroutine remove_negative_rounding_errors_in_slice(y_local_index, x_prev, y_prev, field, local_grid)
Removes the negative rounding errors from a slice of a single field. This works two columns behind an...
Definition: stepfields.F90:433
stepfields_mod::perform_timesmooth_for_field
subroutine perform_timesmooth_for_field(field, zfield, local_grid, x_index, y_index, c1, c2)
Performs initial timesmoothing for a theta or Q field using Robert filter. This is finished off in sw...
Definition: stepfields.F90:516
q_indices_mod::get_q_index
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
stepfields_mod::timestep_callback
subroutine timestep_callback(current_state)
Called at each timestep and will perform swapping and smoothing as required.
Definition: stepfields.F90:237
prognostics_mod::prognostic_field_type
A prognostic field which is assumed to be 3D.
Definition: prognostics.F90:13
stepfields_mod::tend_pr_tot_tabs
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_tabs
Definition: stepfields.F90:38
stepfields_mod::tend_3d_tabs
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_tabs
Definition: stepfields.F90:30
datadefn_mod
Contains common definitions for the data and datatypes used by MONC.
Definition: datadefn.F90:2
stepfields_mod
Does the field stepping Stepping is called at the end of processing a column and steps the x-2 column...
Definition: stepfields.F90:3
stepfields_mod::iqg
integer iqg
Definition: stepfields.F90:46
stepfields_mod::field_information_retrieval_callback
subroutine field_information_retrieval_callback(current_state, name, field_information)
Field information retrieval callback, this returns information for a specific component's published f...
Definition: stepfields.F90:694
q_indices_mod
This manages the Q variables and specifically the mapping between names and the index that they are s...
Definition: q_indices.F90:2
stepfields_mod::l_nonconservative_positive_q
logical l_nonconservative_positive_q
Definition: stepfields.F90:26
stepfields_mod::tend_pr_tot_qs
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_qs
Definition: stepfields.F90:38
registry_mod
MONC component registry.
Definition: registry.F90:5
stepfields_mod::l_tend_pr_tot_th
logical l_tend_pr_tot_th
Definition: stepfields.F90:42
stepfields_mod::l_tend_3d_qg
logical l_tend_3d_qg
Definition: stepfields.F90:34
stepfields_mod::tend_3d_th
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_th
Definition: stepfields.F90:30
stepfields_mod::diagnostic_generation_frequency
integer diagnostic_generation_frequency
Definition: stepfields.F90:47
grids_mod
Functionality to support the different types of grid and abstraction between global grids and local o...
Definition: grids.F90:5
stepfields_mod::iqi
integer iqi
Definition: stepfields.F90:46
optionsdatabase_mod
Manages the options database. Contains administration functions and deduce runtime options from the c...
Definition: optionsdatabase.F90:7
stepfields_mod::initialisation_callback
subroutine initialisation_callback(current_state)
Initialisation callback.
Definition: stepfields.F90:90
monc_component_mod::component_descriptor_type
Description of a component.
Definition: monc_component.F90:42
monc_component_mod::component_double_data_type
integer, parameter, public component_double_data_type
Definition: monc_component.F90:16
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
monc_component_mod::component_array_field_type
integer, parameter, public component_array_field_type
Definition: monc_component.F90:15