MONC
pwadvection.F90
Go to the documentation of this file.
1 
8  use collections_mod, only : map_type
9  use state_mod, only : model_state_type
10  use grids_mod, only : z_index, y_index, x_index
12 
13 implicit none
14 
15 #ifndef TEST_MODE
16  private
17 #endif
18 
20 
21  logical :: l_toplevel=.true.
22 
23  ! Local tendency diagnostic variables for this component
24  ! 3D tendency fields and logicals for their use
25  real(kind=default_precision), dimension(:,:,:), allocatable :: &
32  ! Local mean tendency profile fields and logicals for their use
33  real(kind=default_precision), dimension(:), allocatable :: &
40  ! q indices
41  integer :: iqv=0, iql=0, iqr=0, iqi=0, iqs=0, iqg=0
43 
44 
46 contains
47 
51  pwadvection_get_descriptor%name="pw_advection"
52  pwadvection_get_descriptor%version=0.1
55 
58  allocate(pwadvection_get_descriptor%published_fields(11+11))
59 
60  pwadvection_get_descriptor%published_fields(1)= "tend_u_pwadvection_3d_local"
61  pwadvection_get_descriptor%published_fields(2)= "tend_v_pwadvection_3d_local"
62  pwadvection_get_descriptor%published_fields(3)= "tend_w_pwadvection_3d_local"
63  pwadvection_get_descriptor%published_fields(4)= "tend_th_pwadvection_3d_local"
64  pwadvection_get_descriptor%published_fields(5)= "tend_qv_pwadvection_3d_local"
65  pwadvection_get_descriptor%published_fields(6)= "tend_ql_pwadvection_3d_local"
66  pwadvection_get_descriptor%published_fields(7)= "tend_qi_pwadvection_3d_local"
67  pwadvection_get_descriptor%published_fields(8)= "tend_qr_pwadvection_3d_local"
68  pwadvection_get_descriptor%published_fields(9)= "tend_qs_pwadvection_3d_local"
69  pwadvection_get_descriptor%published_fields(10)="tend_qg_pwadvection_3d_local"
70  pwadvection_get_descriptor%published_fields(11)="tend_tabs_pwadvection_3d_local"
71 
72  pwadvection_get_descriptor%published_fields(11+1)= "tend_u_pwadvection_profile_total_local"
73  pwadvection_get_descriptor%published_fields(11+2)= "tend_v_pwadvection_profile_total_local"
74  pwadvection_get_descriptor%published_fields(11+3)= "tend_w_pwadvection_profile_total_local"
75  pwadvection_get_descriptor%published_fields(11+4)= "tend_th_pwadvection_profile_total_local"
76  pwadvection_get_descriptor%published_fields(11+5)= "tend_qv_pwadvection_profile_total_local"
77  pwadvection_get_descriptor%published_fields(11+6)= "tend_ql_pwadvection_profile_total_local"
78  pwadvection_get_descriptor%published_fields(11+7)= "tend_qi_pwadvection_profile_total_local"
79  pwadvection_get_descriptor%published_fields(11+8)= "tend_qr_pwadvection_profile_total_local"
80  pwadvection_get_descriptor%published_fields(11+9)= "tend_qs_pwadvection_profile_total_local"
81  pwadvection_get_descriptor%published_fields(11+10)="tend_qg_pwadvection_profile_total_local"
82  pwadvection_get_descriptor%published_fields(11+11)="tend_tabs_pwadvection_profile_total_local"
83 
84  end function pwadvection_get_descriptor
85 
88  subroutine initialisation_callback(current_state)
89  type(model_state_type), target, intent(inout) :: current_state
90  logical :: l_qdiag
91 
92  advect_flow=determine_if_advection_here(options_get_string(current_state%options_database, "advection_flow_fields"))
93  advect_th=determine_if_advection_here(options_get_string(current_state%options_database, "advection_theta_field"))
94  advect_q=determine_if_advection_here(options_get_string(current_state%options_database, "advection_q_fields"))
95 
96  ! Set tendency diagnostic logicals based on availability
97  ! Need to use 3d tendencies to compute the profiles, so they will be allocated
98  ! in the case where profiles are available
99  l_qdiag = (.not. current_state%passive_q .and. current_state%number_q_fields .gt. 0) .and. advect_q
100 
101  l_tend_pr_tot_u = current_state%u%active .and. advect_flow
102  l_tend_pr_tot_v = current_state%v%active .and. advect_flow
103  l_tend_pr_tot_w = current_state%w%active .and. advect_flow
104  l_tend_pr_tot_th = current_state%th%active .and. advect_th
105  l_tend_pr_tot_qv = l_qdiag .and. current_state%number_q_fields .ge. 1
106  l_tend_pr_tot_ql = l_qdiag .and. current_state%number_q_fields .ge. 2
107  l_tend_pr_tot_qi = l_qdiag .and. current_state%number_q_fields .ge. 11
108  l_tend_pr_tot_qr = l_qdiag .and. current_state%number_q_fields .ge. 11
109  l_tend_pr_tot_qs = l_qdiag .and. current_state%number_q_fields .ge. 11
110  l_tend_pr_tot_qg = l_qdiag .and. current_state%number_q_fields .ge. 11
112 
113  l_tend_3d_u = (current_state%u%active .and. advect_flow) .or. l_tend_pr_tot_u
114  l_tend_3d_v = (current_state%v%active .and. advect_flow) .or. l_tend_pr_tot_v
115  l_tend_3d_w = (current_state%w%active .and. advect_flow) .or. l_tend_pr_tot_w
116  l_tend_3d_th = (current_state%th%active .and. advect_th) .or. l_tend_pr_tot_th
117  l_tend_3d_qv = (l_qdiag .and. current_state%number_q_fields .ge. 1) .or. l_tend_pr_tot_qv
118  l_tend_3d_ql = (l_qdiag .and. current_state%number_q_fields .ge. 2) .or. l_tend_pr_tot_ql
119  l_tend_3d_qi = (l_qdiag .and. current_state%number_q_fields .ge. 11) .or. l_tend_pr_tot_qi
120  l_tend_3d_qr = (l_qdiag .and. current_state%number_q_fields .ge. 11) .or. l_tend_pr_tot_qr
121  l_tend_3d_qs = (l_qdiag .and. current_state%number_q_fields .ge. 11) .or. l_tend_pr_tot_qs
122  l_tend_3d_qg = (l_qdiag .and. current_state%number_q_fields .ge. 11) .or. l_tend_pr_tot_qg
124 
125  ! Allocate 3d tendency fields upon availability
126  if (l_tend_3d_u) then
127  allocate( tend_3d_u(current_state%local_grid%size(z_index), &
128  current_state%local_grid%size(y_index), &
129  current_state%local_grid%size(x_index) ) )
130  endif
131  if (l_tend_3d_v) then
132  allocate( tend_3d_v(current_state%local_grid%size(z_index), &
133  current_state%local_grid%size(y_index), &
134  current_state%local_grid%size(x_index) ) )
135  endif
136  if (l_tend_3d_w) then
137  allocate( tend_3d_w(current_state%local_grid%size(z_index), &
138  current_state%local_grid%size(y_index), &
139  current_state%local_grid%size(x_index) ) )
140  endif
141  if (l_tend_3d_th) then
142  allocate( tend_3d_th(current_state%local_grid%size(z_index), &
143  current_state%local_grid%size(y_index), &
144  current_state%local_grid%size(x_index) ) )
145  endif
146  if (l_tend_3d_qv) then
147  iqv=get_q_index(standard_q_names%VAPOUR, 'pw_advection')
148  allocate( tend_3d_qv(current_state%local_grid%size(z_index), &
149  current_state%local_grid%size(y_index), &
150  current_state%local_grid%size(x_index) ) )
151  endif
152  if (l_tend_3d_ql) then
153  iql=get_q_index(standard_q_names%CLOUD_LIQUID_MASS, 'pw_advection')
154  allocate( tend_3d_ql(current_state%local_grid%size(z_index), &
155  current_state%local_grid%size(y_index), &
156  current_state%local_grid%size(x_index) ) )
157  endif
158  if (l_tend_3d_qi) then
159  iqi=get_q_index(standard_q_names%ICE_MASS, 'pw_advection')
160  allocate( tend_3d_qi(current_state%local_grid%size(z_index), &
161  current_state%local_grid%size(y_index), &
162  current_state%local_grid%size(x_index) ) )
163  endif
164  if (l_tend_3d_qr) then
165  iqr=get_q_index(standard_q_names%RAIN_MASS, 'pw_advection')
166  allocate( tend_3d_qr(current_state%local_grid%size(z_index), &
167  current_state%local_grid%size(y_index), &
168  current_state%local_grid%size(x_index) ) )
169  endif
170  if (l_tend_3d_qs) then
171  iqs=get_q_index(standard_q_names%SNOW_MASS, 'pw_advection')
172  allocate( tend_3d_qs(current_state%local_grid%size(z_index), &
173  current_state%local_grid%size(y_index), &
174  current_state%local_grid%size(x_index) ) )
175  endif
176  if (l_tend_3d_qg) then
177  iqg=get_q_index(standard_q_names%GRAUPEL_MASS, 'pw_advection')
178  allocate( tend_3d_qg(current_state%local_grid%size(z_index), &
179  current_state%local_grid%size(y_index), &
180  current_state%local_grid%size(x_index) ) )
181  endif
182  if (l_tend_3d_tabs) then
183  allocate( tend_3d_tabs(current_state%local_grid%size(z_index), &
184  current_state%local_grid%size(y_index), &
185  current_state%local_grid%size(x_index) ) )
186  endif
187 
188  ! Allocate profile tendency fields upon availability
189  if (l_tend_pr_tot_u) then
190  allocate( tend_pr_tot_u(current_state%local_grid%size(z_index)) )
191  endif
192  if (l_tend_pr_tot_v) then
193  allocate( tend_pr_tot_v(current_state%local_grid%size(z_index)) )
194  endif
195  if (l_tend_pr_tot_w) then
196  allocate( tend_pr_tot_w(current_state%local_grid%size(z_index)) )
197  endif
198  if (l_tend_pr_tot_th) then
199  allocate( tend_pr_tot_th(current_state%local_grid%size(z_index)) )
200  endif
201  if (l_tend_pr_tot_qv) then
202  allocate( tend_pr_tot_qv(current_state%local_grid%size(z_index)) )
203  endif
204  if (l_tend_pr_tot_ql) then
205  allocate( tend_pr_tot_ql(current_state%local_grid%size(z_index)) )
206  endif
207  if (l_tend_pr_tot_qi) then
208  allocate( tend_pr_tot_qi(current_state%local_grid%size(z_index)) )
209  endif
210  if (l_tend_pr_tot_qr) then
211  allocate( tend_pr_tot_qr(current_state%local_grid%size(z_index)) )
212  endif
213  if (l_tend_pr_tot_qs) then
214  allocate( tend_pr_tot_qs(current_state%local_grid%size(z_index)) )
215  endif
216  if (l_tend_pr_tot_qg) then
217  allocate( tend_pr_tot_qg(current_state%local_grid%size(z_index)) )
218  endif
219  if (l_tend_pr_tot_tabs) then
220  allocate( tend_pr_tot_tabs(current_state%local_grid%size(z_index)) )
221  endif
222 
223  ! Save the sampling_frequency to force diagnostic calculation on select time steps
224  diagnostic_generation_frequency=options_get_integer(current_state%options_database, "sampling_frequency")
225 
226  end subroutine initialisation_callback
227 
228 
229  subroutine finalisation_callback(current_state)
230  type(model_state_type), target, intent(inout) :: current_state
231 
232  if (allocated(tend_3d_u)) deallocate(tend_3d_u)
233  if (allocated(tend_3d_v)) deallocate(tend_3d_v)
234  if (allocated(tend_3d_w)) deallocate(tend_3d_w)
235  if (allocated(tend_3d_th)) deallocate(tend_3d_th)
236  if (allocated(tend_3d_qv)) deallocate(tend_3d_qv)
237  if (allocated(tend_3d_ql)) deallocate(tend_3d_ql)
238  if (allocated(tend_3d_qi)) deallocate(tend_3d_qi)
239  if (allocated(tend_3d_qr)) deallocate(tend_3d_qr)
240  if (allocated(tend_3d_qs)) deallocate(tend_3d_qs)
241  if (allocated(tend_3d_qg)) deallocate(tend_3d_qg)
242  if (allocated(tend_3d_tabs)) deallocate(tend_3d_tabs)
243 
244  if (allocated(tend_pr_tot_u)) deallocate(tend_pr_tot_u)
245  if (allocated(tend_pr_tot_v)) deallocate(tend_pr_tot_v)
246  if (allocated(tend_pr_tot_w)) deallocate(tend_pr_tot_w)
247  if (allocated(tend_pr_tot_th)) deallocate(tend_pr_tot_th)
248  if (allocated(tend_pr_tot_qv)) deallocate(tend_pr_tot_qv)
249  if (allocated(tend_pr_tot_ql)) deallocate(tend_pr_tot_ql)
250  if (allocated(tend_pr_tot_qi)) deallocate(tend_pr_tot_qi)
251  if (allocated(tend_pr_tot_qr)) deallocate(tend_pr_tot_qr)
252  if (allocated(tend_pr_tot_qs)) deallocate(tend_pr_tot_qs)
253  if (allocated(tend_pr_tot_qg)) deallocate(tend_pr_tot_qg)
254  if (allocated(tend_pr_tot_tabs)) deallocate(tend_pr_tot_tabs)
255 
256  end subroutine finalisation_callback
257 
258 
263  subroutine timestep_callback(current_state)
264  type(model_state_type), target, intent(inout) :: current_state
265 
266  integer :: current_x_index, current_y_index, target_x_index, target_y_index
267 
268  current_x_index=current_state%column_local_x
269  current_y_index=current_state%column_local_y
270  target_y_index=current_y_index-current_state%local_grid%halo_size(y_index)
271  target_x_index=current_x_index-current_state%local_grid%halo_size(x_index)
272 
273 
274  ! Zero profile tendency totals on first instance in the sum
275  if (current_state%first_timestep_column) then
276  if (l_tend_pr_tot_u) then
277  tend_pr_tot_u(:)= 0.0_default_precision
278  endif
279  if (l_tend_pr_tot_v) then
280  tend_pr_tot_v(:)= 0.0_default_precision
281  endif
282  if (l_tend_pr_tot_w) then
283  tend_pr_tot_w(:)= 0.0_default_precision
284  endif
285  if (l_tend_pr_tot_th) then
286  tend_pr_tot_th(:)=0.0_default_precision
287  endif
288  if (l_tend_pr_tot_qv) then
289  tend_pr_tot_qv(:)=0.0_default_precision
290  endif
291  if (l_tend_pr_tot_ql) then
292  tend_pr_tot_ql(:)=0.0_default_precision
293  endif
294  if (l_tend_pr_tot_qi) then
295  tend_pr_tot_qi(:)=0.0_default_precision
296  endif
297  if (l_tend_pr_tot_qr) then
298  tend_pr_tot_qr(:)=0.0_default_precision
299  endif
300  if (l_tend_pr_tot_qs) then
301  tend_pr_tot_qs(:)=0.0_default_precision
302  endif
303  if (l_tend_pr_tot_qg) then
304  tend_pr_tot_qg(:)=0.0_default_precision
305  endif
306  if (l_tend_pr_tot_tabs) then
307  tend_pr_tot_tabs(:)=0.0_default_precision
308  endif
309  endif ! zero totals
310 
311 
312  if (current_state%halo_column) return
313 
314  if (mod(current_state%timestep, diagnostic_generation_frequency) == 0) then
315  call save_precomponent_tendencies(current_state, current_x_index, current_y_index, target_x_index, target_y_index)
316  end if
317 
318  if (advect_flow) call advect_flow_fields(current_state, current_x_index, current_y_index)
319  if (advect_th) call advect_th_field(current_state, current_x_index, current_y_index)
320  if (advect_q) call advect_q_field(current_state, current_x_index, current_y_index)
321 
322  if (mod(current_state%timestep, diagnostic_generation_frequency) == 0) then
323  call compute_component_tendencies(current_state, current_x_index, current_y_index, target_x_index, target_y_index)
324  end if
325 
326  end subroutine timestep_callback
327 
332  subroutine advect_q_field(current_state, current_x_index, current_y_index)
333  type(model_state_type), target, intent(inout) :: current_state
334  integer, intent(in) :: current_x_index, current_y_index
335 
336  integer :: k, n
337 
338  do n=1,current_state%number_q_fields
339  do k=2,current_state%local_grid%size(z_index)-1
340 #ifdef U_ACTIVE
341  current_state%sq(n)%data(k, current_y_index, current_x_index)=&
342  current_state%sq(n)%data(k, current_y_index, current_x_index)+current_state%global_grid%configuration%horizontal%cx*&
343  0.5_default_precision*(current_state%u%data(k, current_y_index, current_x_index-1)*&
344  current_state%q(n)%data(k, current_y_index, current_x_index-1)-&
345  current_state%u%data(k, current_y_index, current_x_index)*&
346  current_state%q(n)%data(k, current_y_index, current_x_index+1))
347 #endif
348 #ifdef V_ACTIVE
349  current_state%sq(n)%data(k, current_y_index, current_x_index)=&
350  current_state%sq(n)%data(k, current_y_index, current_x_index)+&
351  current_state%global_grid%configuration%horizontal%cy*0.5_default_precision*&
352  (current_state%v%data(k, current_y_index-1, current_x_index)*&
353  current_state%q(n)%data(k, current_y_index-1, current_x_index)-&
354  current_state%v%data(k, current_y_index, current_x_index)*&
355  current_state%q(n)%data(k, current_y_index+1, current_x_index))
356 #endif
357 #ifdef W_ACTIVE
358  current_state%sq(n)%data(k, current_y_index, current_x_index)=&
359  current_state%sq(n)%data(k, current_y_index, current_x_index)+&
360  2.0_default_precision*(current_state%global_grid%configuration%vertical%tzc1(k)*&
361  current_state%w%data(k-1, current_y_index, current_x_index)*&
362  current_state%q(n)%data(k-1, current_y_index, current_x_index)-&
363  current_state%global_grid%configuration%vertical%tzc2(k)*&
364  current_state%w%data(k, current_y_index, current_x_index)*&
365  current_state%q(n)%data(k+1, current_y_index, current_x_index))
366 #endif
367  end do
368 
369  if (l_toplevel)then
370  k=current_state%local_grid%size(z_index)
371 #ifdef U_ACTIVE
372  current_state%sq(n)%data(k, current_y_index, current_x_index)=&
373  current_state%sq(n)%data(k, current_y_index, current_x_index)+current_state%global_grid%configuration%horizontal%cx*&
374  0.5_default_precision*(current_state%u%data(k, current_y_index, current_x_index-1)*&
375  current_state%q(n)%data(k, current_y_index, current_x_index-1)-&
376  current_state%u%data(k, current_y_index, current_x_index)*&
377  current_state%q(n)%data(k, current_y_index, current_x_index+1))
378 #endif
379 #ifdef V_ACTIVE
380  current_state%sq(n)%data(k, current_y_index, current_x_index)=&
381  current_state%sq(n)%data(k, current_y_index, current_x_index)+current_state%global_grid%configuration%horizontal%cy*&
382  0.5_default_precision*(current_state%v%data(k, current_y_index-1, current_x_index)*&
383  current_state%q(n)%data(k, current_y_index-1, current_x_index)-&
384  current_state%v%data(k, current_y_index, current_x_index)*&
385  current_state%q(n)%data(k, current_y_index+1, current_x_index))
386 #endif
387 #ifdef W_ACTIVE
388  current_state%sq(n)%data(k, current_y_index, current_x_index)=&
389  current_state%sq(n)%data(k, current_y_index, current_x_index)+&
390  current_state%global_grid%configuration%vertical%tzc1(k)*2.0_default_precision*&
391  current_state%w%data(k-1, current_y_index, current_x_index)*&
392  current_state%q(n)%data(k-1, current_y_index, current_x_index)
393 #endif
394  end if
395  end do
396  end subroutine advect_q_field
397 
402  subroutine advect_th_field(current_state, current_x_index, current_y_index)
403  type(model_state_type), target, intent(inout) :: current_state
404  integer, intent(in) :: current_x_index, current_y_index
405 
406  integer :: k
407 
408  if (current_state%th%active) then
409 
410  do k=2,current_state%local_grid%size(z_index)-1
411 #ifdef U_ACTIVE
412  current_state%sth%data(k, current_y_index, current_x_index)= & !current_state%sth%data(k, current_y_index, current_x_index)+&
413  current_state%global_grid%configuration%horizontal%cx*&
414  0.5_default_precision*(current_state%u%data(k, current_y_index, current_x_index-1)*&
415  current_state%th%data(k, current_y_index, current_x_index-1)-&
416  current_state%u%data(k, current_y_index, current_x_index)*&
417  current_state%th%data(k, current_y_index, current_x_index+1))
418 #endif
419 #ifdef V_ACTIVE
420  current_state%sth%data(k, current_y_index, current_x_index)=current_state%sth%data(k, current_y_index, current_x_index)+&
421  current_state%global_grid%configuration%horizontal%cy*0.5_default_precision*&
422  (current_state%v%data(k, current_y_index-1, current_x_index)*&
423  current_state%th%data(k, current_y_index-1, current_x_index)-&
424  current_state%v%data(k, current_y_index, current_x_index)*&
425  current_state%th%data(k, current_y_index+1, current_x_index))
426 #endif
427 #ifdef W_ACTIVE
428  current_state%sth%data(k, current_y_index, current_x_index)=current_state%sth%data(k, current_y_index, current_x_index)+&
429  2.0_default_precision*(current_state%global_grid%configuration%vertical%tzc1(k)*&
430  current_state%w%data(k-1, current_y_index, current_x_index)*&
431  current_state%th%data(k-1, current_y_index, current_x_index)-&
432  current_state%global_grid%configuration%vertical%tzc2(k)*&
433  current_state%w%data(k, current_y_index, current_x_index)*&
434  current_state%th%data(k+1, current_y_index, current_x_index))
435 #endif
436  end do
437 
438  if (l_toplevel)then
439 
440  k=current_state%local_grid%size(z_index)
441 #ifdef U_ACTIVE
442  current_state%sth%data(k, current_y_index, current_x_index)= & !current_state%sth%data(k, current_y_index, current_x_index)+&
443  current_state%global_grid%configuration%horizontal%cx*&
444  0.5_default_precision*(current_state%u%data(k, current_y_index, current_x_index-1)*&
445  current_state%th%data(k, current_y_index, current_x_index-1)-&
446  current_state%u%data(k, current_y_index, current_x_index)*&
447  current_state%th%data(k, current_y_index, current_x_index+1))
448 #endif
449 #ifdef V_ACTIVE
450  current_state%sth%data(k, current_y_index, current_x_index)=current_state%sth%data(k, current_y_index, current_x_index)+&
451  current_state%global_grid%configuration%horizontal%cy*0.5_default_precision*&
452  (current_state%v%data(k, current_y_index-1, current_x_index)*&
453  current_state%th%data(k, current_y_index-1, current_x_index)-&
454  current_state%v%data(k, current_y_index, current_x_index)*&
455  current_state%th%data(k, current_y_index+1, current_x_index))
456 #endif
457 #ifdef W_ACTIVE
458  current_state%sth%data(k, current_y_index, current_x_index)=current_state%sth%data(k, current_y_index, current_x_index)+&
459  current_state%global_grid%configuration%vertical%tzc1(k)*2.0_default_precision*&
460  current_state%w%data(k-1, current_y_index, current_x_index)*current_state%th%data(k-1, current_y_index, &
461  current_x_index)
462 #endif
463  end if
464  end if
465  end subroutine advect_th_field
466 
471  subroutine advect_flow_fields(current_state, current_x_index, current_y_index)
472  type(model_state_type), target, intent(inout) :: current_state
473  integer, intent(in) :: current_x_index, current_y_index
474 
475  integer :: k
476 
477  do k=2,current_state%local_grid%size(z_index)-1
478 #ifdef U_ACTIVE
479  current_state%su%data(k, current_y_index, current_x_index)=&
480  current_state%global_grid%configuration%horizontal%tcx*(current_state%u%data(k, current_y_index, current_x_index-1)*&
481  (current_state%u%data(k, current_y_index, current_x_index)+&
482  current_state%u%data(k, current_y_index, current_x_index-1))-&
483  current_state%u%data(k, current_y_index, current_x_index+1)*&
484  (current_state%u%data(k, current_y_index, current_x_index)+&
485  current_state%u%data(k, current_y_index, current_x_index+1)))
486 #ifdef V_ACTIVE
487  current_state%su%data(k, current_y_index, current_x_index)=current_state%su%data(k, current_y_index, current_x_index)+&
488  current_state%global_grid%configuration%horizontal%tcy*(current_state%u%data(k, current_y_index-1, current_x_index)*&
489  (current_state%v%data(k, current_y_index-1, current_x_index)+&
490  current_state%v%data(k, current_y_index-1, current_x_index+1))-&
491  current_state%u%data(k, current_y_index+1, current_x_index)*&
492  (current_state%v%data(k, current_y_index, current_x_index)+&
493  current_state%v%data(k, current_y_index, current_x_index+1)))
494 #endif
495 #ifdef W_ACTIVE
496  current_state%su%data(k, current_y_index, current_x_index)=current_state%su%data(k, current_y_index, current_x_index)+&
497  (current_state%global_grid%configuration%vertical%tzc1(k)*current_state%u%data(k-1, current_y_index, current_x_index)*&
498  (current_state%w%data(k-1, current_y_index, current_x_index)+&
499  current_state%w%data(k-1, current_y_index, current_x_index+1))-&
500  current_state%global_grid%configuration%vertical%tzc2(k)*current_state%u%data(k+1, current_y_index, current_x_index)*&
501  (current_state%w%data(k, current_y_index, current_x_index)+&
502  current_state%w%data(k, current_y_index, current_x_index+1)))
503 #endif
504 #endif
505 
506 #ifdef V_ACTIVE
507  current_state%sv%data(k, current_y_index, current_x_index)=current_state%global_grid%configuration%horizontal%tcy*(&
508  current_state%v%data(k, current_y_index-1, current_x_index)*&
509  (current_state%v%data(k, current_y_index, current_x_index)+&
510  current_state%v%data(k, current_y_index-1, current_x_index))-&
511  current_state%v%data(k, current_y_index+1, current_x_index)*&
512  (current_state%v%data(k, current_y_index, current_x_index)+&
513  current_state%v%data(k, current_y_index+1, current_x_index)))
514 #ifdef U_ACTIVE
515  current_state%sv%data(k, current_y_index, current_x_index)=current_state%sv%data(k, current_y_index, current_x_index)+&
516  current_state%global_grid%configuration%horizontal%tcx*(current_state%v%data(k, current_y_index, current_x_index-1)*&
517  (current_state%u%data(k, current_y_index, current_x_index-1)+&
518  current_state%u%data(k, current_y_index+1, current_x_index-1))-&
519  current_state%v%data(k, current_y_index, current_x_index+1)*&
520  (current_state%u%data(k, current_y_index, current_x_index)+&
521  current_state%u%data(k, current_y_index+1, current_x_index)))
522 #endif
523 #ifdef W_ACTIVE
524  current_state%sv%data(k, current_y_index, current_x_index)=current_state%sv%data(k, current_y_index, current_x_index)+&
525  (current_state%global_grid%configuration%vertical%tzc1(k)*current_state%v%data(k-1, current_y_index, current_x_index)*&
526  (current_state%w%data(k-1, current_y_index, current_x_index)+&
527  current_state%w%data(k-1, current_y_index+1, current_x_index))-&
528  current_state%global_grid%configuration%vertical%tzc2(k)*current_state%v%data(k+1, current_y_index, current_x_index)*&
529  (current_state%w%data(k, current_y_index, current_x_index)+&
530  current_state%w%data(k, current_y_index+1, current_x_index)))
531 #endif
532 #endif
533 
534 #ifdef W_ACTIVE
535  current_state%sw%data(k, current_y_index, current_x_index)=(current_state%global_grid%configuration%vertical%tzd1(k)*&
536  current_state%w%data(k-1, current_y_index, current_x_index)*&
537  (current_state%w%data(k, current_y_index, current_x_index)+&
538  current_state%w%data(k-1, current_y_index, current_x_index))-&
539  current_state%global_grid%configuration%vertical%tzd2(k)*current_state%w%data(k+1, current_y_index, current_x_index)*&
540  (current_state%w%data(k, current_y_index, current_x_index)+&
541  current_state%w%data(k+1, current_y_index, current_x_index)))
542 #ifdef U_ACTIVE
543  current_state%sw%data(k, current_y_index, current_x_index)=current_state%sw%data(k, current_y_index, current_x_index)+&
544  current_state%global_grid%configuration%horizontal%tcx*(current_state%w%data(k, current_y_index, current_x_index-1)*&
545  (current_state%u%data(k, current_y_index, current_x_index-1)+&
546  current_state%u%data(k+1, current_y_index, current_x_index-1))-&
547  current_state%w%data(k, current_y_index, current_x_index+1)*&
548  (current_state%u%data(k, current_y_index, current_x_index)+&
549  current_state%u%data(k+1, current_y_index, current_x_index)))
550 #endif
551 #ifdef V_ACTIVE
552  current_state%sw%data(k, current_y_index, current_x_index)=current_state%sw%data(k, current_y_index, current_x_index)+&
553  current_state%global_grid%configuration%horizontal%tcy*(current_state%w%data(k, current_y_index-1, current_x_index)*&
554  (current_state%v%data(k, current_y_index-1, current_x_index)+&
555  current_state%v%data(k+1, current_y_index-1, current_x_index))-&
556  current_state%w%data(k, current_y_index+1, current_x_index)*&
557  (current_state%v%data(k, current_y_index, current_x_index)+&
558  current_state%v%data(k+1, current_y_index, current_x_index)))
559 #endif
560 #endif
561  end do
562 
563  if (l_toplevel)then
564  k=current_state%local_grid%size(z_index)
565 #ifdef U_ACTIVE
566  current_state%su%data(k, current_y_index, current_x_index)=current_state%global_grid%configuration%horizontal%tcx*&
567  (current_state%u%data(k, current_y_index, current_x_index-1)*&
568  (current_state%u%data(k, current_y_index, current_x_index)+&
569  current_state%u%data(k, current_y_index, current_x_index-1))-&
570  current_state%u%data(k, current_y_index, current_x_index+1)*&
571  (current_state%u%data(k, current_y_index, current_x_index)+&
572  current_state%u%data(k, current_y_index, current_x_index+1)))
573 #ifdef V_ACTIVE
574  current_state%su%data(k, current_y_index, current_x_index)=current_state%su%data(k, current_y_index, current_x_index)+&
575  current_state%global_grid%configuration%horizontal%tcy*(current_state%u%data(k, current_y_index-1, current_x_index)*&
576  (current_state%v%data(k, current_y_index-1, current_x_index)+&
577  current_state%v%data(k, current_y_index-1, current_x_index+1))-&
578  current_state%u%data(k, current_y_index+1, current_x_index)*&
579  (current_state%v%data(k, current_y_index, current_x_index)+&
580  current_state%v%data(k, current_y_index, current_x_index+1)))
581 #endif
582 #ifdef W_ACTIVE
583  current_state%su%data(k, current_y_index, current_x_index)=current_state%su%data(k, current_y_index, current_x_index)+&
584  current_state%global_grid%configuration%vertical%tzc1(k)*current_state%u%data(k-1, current_y_index, current_x_index)*&
585  (current_state%w%data(k-1, current_y_index, current_x_index)+&
586  current_state%w%data(k-1, current_y_index, current_x_index+1))
587 #endif
588 #endif
589 
590 #ifdef V_ACTIVE
591  current_state%sv%data(k, current_y_index, current_x_index)=current_state%global_grid%configuration%horizontal%tcy*&
592  (current_state%v%data(k, current_y_index-1, current_x_index)*&
593  (current_state%v%data(k, current_y_index, current_x_index)+&
594  current_state%v%data(k, current_y_index-1, current_x_index))-&
595  current_state%v%data(k, current_y_index+1, current_x_index)*&
596  (current_state%v%data(k, current_y_index, current_x_index)+&
597  current_state%v%data(k, current_y_index+1, current_x_index)))
598 #ifdef U_ACTIVE
599  current_state%sv%data(k, current_y_index, current_x_index)=current_state%sv%data(k, current_y_index, current_x_index)+&
600  current_state%global_grid%configuration%horizontal%tcx*(current_state%v%data(k, current_y_index, current_x_index-1)*&
601  (current_state%u%data(k, current_y_index, current_x_index-1)+&
602  current_state%u%data(k, current_y_index+1, current_x_index-1))-&
603  current_state%v%data(k, current_y_index, current_x_index+1)*&
604  (current_state%u%data(k, current_y_index, current_x_index)+&
605  current_state%u%data(k, current_y_index+1, current_x_index)))
606 #endif
607 #ifdef W_ACTIVE
608  current_state%sv%data(k, current_y_index, current_x_index)=current_state%sv%data(k, current_y_index, current_x_index)+&
609  current_state%global_grid%configuration%vertical%tzc1(k)*current_state%v%data(k-1, current_y_index, current_x_index)*&
610  (current_state%w%data(k-1, current_y_index, current_x_index)+&
611  current_state%w%data(k-1, current_y_index+1, current_x_index))
612 #endif
613 #endif
614  end if
615  end subroutine advect_flow_fields
616 
617 
624  subroutine save_precomponent_tendencies(current_state, cxn, cyn, txn, tyn)
625  type(model_state_type), target, intent(in) :: current_state
626  integer, intent(in) :: cxn, cyn, txn, tyn
627 
628  ! Save 3d tendency fields upon request (of 3d or profiles) and availability
629  if (l_tend_3d_u) then
630  tend_3d_u(:,tyn,txn)=current_state%su%data(:,cyn,cxn)
631  endif
632  if (l_tend_3d_v) then
633  tend_3d_v(:,tyn,txn)=current_state%sv%data(:,cyn,cxn)
634  endif
635  if (l_tend_3d_w) then
636  tend_3d_w(:,tyn,txn)=current_state%sw%data(:,cyn,cxn)
637  endif
638  if (l_tend_3d_th) then
639  tend_3d_th(:,tyn,txn)=current_state%sth%data(:,cyn,cxn)
640  endif
641  if (l_tend_3d_qv) then
642  tend_3d_qv(:,tyn,txn)=current_state%sq(iqv)%data(:,cyn,cxn)
643  endif
644  if (l_tend_3d_ql) then
645  tend_3d_ql(:,tyn,txn)=current_state%sq(iql)%data(:,cyn,cxn)
646  endif
647  if (l_tend_3d_qi) then
648  tend_3d_qi(:,tyn,txn)=current_state%sq(iqi)%data(:,cyn,cxn)
649  endif
650  if (l_tend_3d_qr) then
651  tend_3d_qr(:,tyn,txn)=current_state%sq(iqr)%data(:,cyn,cxn)
652  endif
653  if (l_tend_3d_qs) then
654  tend_3d_qs(:,tyn,txn)=current_state%sq(iqs)%data(:,cyn,cxn)
655  endif
656  if (l_tend_3d_qg) then
657  tend_3d_qg(:,tyn,txn)=current_state%sq(iqg)%data(:,cyn,cxn)
658  endif
659  if (l_tend_3d_tabs) then
660  tend_3d_tabs(:,tyn,txn)=current_state%sth%data(:,cyn,cxn) * current_state%global_grid%configuration%vertical%rprefrcp(:)
661  endif
662 
663  end subroutine save_precomponent_tendencies
664 
665 
672  subroutine compute_component_tendencies(current_state, cxn, cyn, txn, tyn)
673  type(model_state_type), target, intent(inout) :: current_state
674  integer, intent(in) :: cxn, cyn, txn, tyn
675 
676  ! Calculate change in tendency due to component
677  if (l_tend_3d_u) then
678  tend_3d_u(:,tyn,txn)=current_state%su%data(:,cyn,cxn) - tend_3d_u(:,tyn,txn)
679  endif
680  if (l_tend_3d_v) then
681  tend_3d_v(:,tyn,txn)=current_state%sv%data(:,cyn,cxn) - tend_3d_v(:,tyn,txn)
682  endif
683  if (l_tend_3d_w) then
684  tend_3d_w(:,tyn,txn)=current_state%sw%data(:,cyn,cxn) - tend_3d_w(:,tyn,txn)
685  endif
686  if (l_tend_3d_th) then
687  tend_3d_th(:,tyn,txn)=current_state%sth%data(:,cyn,cxn) - tend_3d_th(:,tyn,txn)
688  endif
689  if (l_tend_3d_qv) then
690  tend_3d_qv(:,tyn,txn)=current_state%sq(iqv)%data(:,cyn,cxn) - tend_3d_qv(:,tyn,txn)
691  endif
692  if (l_tend_3d_ql) then
693  tend_3d_ql(:,tyn,txn)=current_state%sq(iql)%data(:,cyn,cxn) - tend_3d_ql(:,tyn,txn)
694  endif
695  if (l_tend_3d_qi) then
696  tend_3d_qi(:,tyn,txn)=current_state%sq(iqi)%data(:,cyn,cxn) - tend_3d_qi(:,tyn,txn)
697  endif
698  if (l_tend_3d_qr) then
699  tend_3d_qr(:,tyn,txn)=current_state%sq(iqr)%data(:,cyn,cxn) - tend_3d_qr(:,tyn,txn)
700  endif
701  if (l_tend_3d_qs) then
702  tend_3d_qs(:,tyn,txn)=current_state%sq(iqs)%data(:,cyn,cxn) - tend_3d_qs(:,tyn,txn)
703  endif
704  if (l_tend_3d_qg) then
705  tend_3d_qg(:,tyn,txn)=current_state%sq(iqg)%data(:,cyn,cxn) - tend_3d_qg(:,tyn,txn)
706  endif
707  if (l_tend_3d_tabs) then
708  tend_3d_tabs(:,tyn,txn)= &
709  current_state%sth%data(:,cyn,cxn) * current_state%global_grid%configuration%vertical%rprefrcp(:) &
710  - tend_3d_tabs(:,tyn,txn)
711  endif
712 
713  ! Add local tendency fields to the profile total
714  if (l_tend_pr_tot_u) then
715  tend_pr_tot_u(:)=tend_pr_tot_u(:) + tend_3d_u(:,tyn,txn)
716  endif
717  if (l_tend_pr_tot_v) then
718  tend_pr_tot_v(:)=tend_pr_tot_v(:) + tend_3d_v(:,tyn,txn)
719  endif
720  if (l_tend_pr_tot_w) then
721  tend_pr_tot_w(:)=tend_pr_tot_w(:) + tend_3d_w(:,tyn,txn)
722  endif
723  if (l_tend_pr_tot_th) then
724  tend_pr_tot_th(:)=tend_pr_tot_th(:) + tend_3d_th(:,tyn,txn)
725  endif
726  if (l_tend_pr_tot_qv) then
727  tend_pr_tot_qv(:)=tend_pr_tot_qv(:) + tend_3d_qv(:,tyn,txn)
728  endif
729  if (l_tend_pr_tot_ql) then
730  tend_pr_tot_ql(:)=tend_pr_tot_ql(:) + tend_3d_ql(:,tyn,txn)
731  endif
732  if (l_tend_pr_tot_qi) then
733  tend_pr_tot_qi(:)=tend_pr_tot_qi(:) + tend_3d_qi(:,tyn,txn)
734  endif
735  if (l_tend_pr_tot_qr) then
736  tend_pr_tot_qr(:)=tend_pr_tot_qr(:) + tend_3d_qr(:,tyn,txn)
737  endif
738  if (l_tend_pr_tot_qs) then
739  tend_pr_tot_qs(:)=tend_pr_tot_qs(:) + tend_3d_qs(:,tyn,txn)
740  endif
741  if (l_tend_pr_tot_qg) then
742  tend_pr_tot_qg(:)=tend_pr_tot_qg(:) + tend_3d_qg(:,tyn,txn)
743  endif
744  if (l_tend_pr_tot_tabs) then
746  endif
747 
748  end subroutine compute_component_tendencies
749 
750 
756  subroutine field_information_retrieval_callback(current_state, name, field_information)
757  type(model_state_type), target, intent(inout) :: current_state
758  character(len=*), intent(in) :: name
759  type(component_field_information_type), intent(out) :: field_information
760  integer :: strcomp
761 
762  ! Field information for 3d
763  strcomp=index(name, "pwadvection_3d_local")
764  if (strcomp .ne. 0) then
765  field_information%field_type=component_array_field_type
766  field_information%number_dimensions=3
767  field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
768  field_information%dimension_sizes(2)=current_state%local_grid%size(y_index)
769  field_information%dimension_sizes(3)=current_state%local_grid%size(x_index)
770  field_information%data_type=component_double_data_type
771 
772  if (name .eq. "tend_u_pwadvection_3d_local") then
773  field_information%enabled=l_tend_3d_u
774  else if (name .eq. "tend_v_pwadvection_3d_local") then
775  field_information%enabled=l_tend_3d_v
776  else if (name .eq. "tend_w_pwadvection_3d_local") then
777  field_information%enabled=l_tend_3d_w
778  else if (name .eq. "tend_th_pwadvection_3d_local") then
779  field_information%enabled=l_tend_3d_th
780  else if (name .eq. "tend_qv_pwadvection_3d_local") then
781  field_information%enabled=l_tend_3d_qv
782  else if (name .eq. "tend_ql_pwadvection_3d_local") then
783  field_information%enabled=l_tend_3d_ql
784  else if (name .eq. "tend_qi_pwadvection_3d_local") then
785  field_information%enabled=l_tend_3d_qi
786  else if (name .eq. "tend_qr_pwadvection_3d_local") then
787  field_information%enabled=l_tend_3d_qr
788  else if (name .eq. "tend_qs_pwadvection_3d_local") then
789  field_information%enabled=l_tend_3d_qs
790  else if (name .eq. "tend_qg_pwadvection_3d_local") then
791  field_information%enabled=l_tend_3d_qg
792  else if (name .eq. "tend_tabs_pwadvection_3d_local") then
793  field_information%enabled=l_tend_3d_tabs
794  else
795  field_information%enabled=.true.
796  end if
797 
798  end if !end 3d check
799 
800  ! Field information for profiles
801  strcomp=index(name, "pwadvection_profile_total_local")
802  if (strcomp .ne. 0) then
803  field_information%field_type=component_array_field_type
804  field_information%number_dimensions=1
805  field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
806  field_information%data_type=component_double_data_type
807 
808  if (name .eq. "tend_u_pwadvection_profile_total_local") then
809  field_information%enabled=l_tend_pr_tot_u
810  else if (name .eq. "tend_v_pwadvection_profile_total_local") then
811  field_information%enabled=l_tend_pr_tot_v
812  else if (name .eq. "tend_w_pwadvection_profile_total_local") then
813  field_information%enabled=l_tend_pr_tot_w
814  else if (name .eq. "tend_th_pwadvection_profile_total_local") then
815  field_information%enabled=l_tend_pr_tot_th
816  else if (name .eq. "tend_qv_pwadvection_profile_total_local") then
817  field_information%enabled=l_tend_pr_tot_qv
818  else if (name .eq. "tend_ql_pwadvection_profile_total_local") then
819  field_information%enabled=l_tend_pr_tot_ql
820  else if (name .eq. "tend_qi_pwadvection_profile_total_local") then
821  field_information%enabled=l_tend_pr_tot_qi
822  else if (name .eq. "tend_qr_pwadvection_profile_total_local") then
823  field_information%enabled=l_tend_pr_tot_qr
824  else if (name .eq. "tend_qs_pwadvection_profile_total_local") then
825  field_information%enabled=l_tend_pr_tot_qs
826  else if (name .eq. "tend_qg_pwadvection_profile_total_local") then
827  field_information%enabled=l_tend_pr_tot_qg
828  else if (name .eq. "tend_tabs_pwadvection_profile_total_local") then
829  field_information%enabled=l_tend_pr_tot_tabs
830  else
831  field_information%enabled=.true.
832  end if
833 
834  end if !end profile check
835 
837 
838 
843  subroutine field_value_retrieval_callback(current_state, name, field_value)
844  type(model_state_type), target, intent(inout) :: current_state
845  character(len=*), intent(in) :: name
846  type(component_field_value_type), intent(out) :: field_value
847 
848  ! 3d Tendency Fields
849  if (name .eq. "tend_u_pwadvection_3d_local" .and. allocated(tend_3d_u)) then
850  call set_published_field_value(field_value, real_3d_field=tend_3d_u)
851  else if (name .eq. "tend_v_pwadvection_3d_local" .and. allocated(tend_3d_v)) then
852  call set_published_field_value(field_value, real_3d_field=tend_3d_v)
853  else if (name .eq. "tend_w_pwadvection_3d_local" .and. allocated(tend_3d_w)) then
854  call set_published_field_value(field_value, real_3d_field=tend_3d_w)
855  else if (name .eq. "tend_th_pwadvection_3d_local" .and. allocated(tend_3d_th)) then
856  call set_published_field_value(field_value, real_3d_field=tend_3d_th)
857  else if (name .eq. "tend_qv_pwadvection_3d_local" .and. allocated(tend_3d_qv)) then
858  call set_published_field_value(field_value, real_3d_field=tend_3d_qv)
859  else if (name .eq. "tend_ql_pwadvection_3d_local" .and. allocated(tend_3d_ql)) then
860  call set_published_field_value(field_value, real_3d_field=tend_3d_ql)
861  else if (name .eq. "tend_qi_pwadvection_3d_local" .and. allocated(tend_3d_qi)) then
862  call set_published_field_value(field_value, real_3d_field=tend_3d_qi)
863  else if (name .eq. "tend_qr_pwadvection_3d_local" .and. allocated(tend_3d_qr)) then
864  call set_published_field_value(field_value, real_3d_field=tend_3d_qr)
865  else if (name .eq. "tend_qs_pwadvection_3d_local" .and. allocated(tend_3d_qs)) then
866  call set_published_field_value(field_value, real_3d_field=tend_3d_qs)
867  else if (name .eq. "tend_qg_pwadvection_3d_local" .and. allocated(tend_3d_qg)) then
868  call set_published_field_value(field_value, real_3d_field=tend_3d_qg)
869  else if (name .eq. "tend_tabs_pwadvection_3d_local" .and. allocated(tend_3d_tabs)) then
870  call set_published_field_value(field_value, real_3d_field=tend_3d_tabs)
871 
872  ! Profile Tendency Fields
873  else if (name .eq. "tend_u_pwadvection_profile_total_local" .and. allocated(tend_pr_tot_u)) then
874  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_u)
875  else if (name .eq. "tend_v_pwadvection_profile_total_local" .and. allocated(tend_pr_tot_v)) then
876  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_v)
877  else if (name .eq. "tend_w_pwadvection_profile_total_local" .and. allocated(tend_pr_tot_w)) then
878  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_w)
879  else if (name .eq. "tend_th_pwadvection_profile_total_local" .and. allocated(tend_pr_tot_th)) then
880  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_th)
881  else if (name .eq. "tend_qv_pwadvection_profile_total_local" .and. allocated(tend_pr_tot_qv)) then
882  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_qv)
883  else if (name .eq. "tend_ql_pwadvection_profile_total_local" .and. allocated(tend_pr_tot_ql)) then
884  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_ql)
885  else if (name .eq. "tend_qi_pwadvection_profile_total_local" .and. allocated(tend_pr_tot_qi)) then
886  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_qi)
887  else if (name .eq. "tend_qr_pwadvection_profile_total_local" .and. allocated(tend_pr_tot_qr)) then
888  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_qr)
889  else if (name .eq. "tend_qs_pwadvection_profile_total_local" .and. allocated(tend_pr_tot_qs)) then
890  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_qs)
891  else if (name .eq. "tend_qg_pwadvection_profile_total_local" .and. allocated(tend_pr_tot_qg)) then
892  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_qg)
893  else if (name .eq. "tend_tabs_pwadvection_profile_total_local" .and. allocated(tend_pr_tot_tabs)) then
894  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_tabs)
895  end if
896 
897  end subroutine field_value_retrieval_callback
898 
899 
904  subroutine set_published_field_value(field_value, real_1d_field, real_2d_field, real_3d_field)
905  type(component_field_value_type), intent(inout) :: field_value
906  real(kind=default_precision), dimension(:), optional :: real_1d_field
907  real(kind=default_precision), dimension(:,:), optional :: real_2d_field
908  real(kind=default_precision), dimension(:,:,:), optional :: real_3d_field
909 
910  if (present(real_1d_field)) then
911  allocate(field_value%real_1d_array(size(real_1d_field)), source=real_1d_field)
912  else if (present(real_2d_field)) then
913  allocate(field_value%real_2d_array(size(real_2d_field, 1), size(real_2d_field, 2)), source=real_2d_field)
914  else if (present(real_3d_field)) then
915  allocate(field_value%real_3d_array(size(real_3d_field, 1), size(real_3d_field, 2), size(real_3d_field, 3)), &
916  source=real_3d_field)
917  end if
918  end subroutine set_published_field_value
919 
920 
925  logical function determine_if_advection_here(field)
926  character(len=*), intent(in) :: field
927 
928  if (len_trim(field) .ne. 0) then
929  if (trim(field) .eq. "pw" .or. trim(field) .eq. "any") then
931  else
933  end if
934  else
936  end if
937  end function determine_if_advection_here
938 
939 end module pwadvection_mod
940 
pwadvection_mod::iqg
integer iqg
Definition: pwadvection.F90:41
collections_mod::map_type
Map data structure that holds string (length 20 maximum) key value pairs.
Definition: collections.F90:86
pwadvection_mod
Piacsek-Williams advection scheme.
Definition: pwadvection.F90:2
pwadvection_mod::l_tend_pr_tot_th
logical l_tend_pr_tot_th
Definition: pwadvection.F90:37
pwadvection_mod::l_tend_3d_ql
logical l_tend_3d_ql
Definition: pwadvection.F90:29
pwadvection_mod::tend_3d_qg
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_qg
Definition: pwadvection.F90:25
pwadvection_mod::tend_pr_tot_v
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_v
Definition: pwadvection.F90:33
pwadvection_mod::advect_flow_fields
subroutine advect_flow_fields(current_state, current_x_index, current_y_index)
Advects the flow fields depending upon which fields are active in the model in a column.
Definition: pwadvection.F90:472
pwadvection_mod::l_tend_pr_tot_qv
logical l_tend_pr_tot_qv
Definition: pwadvection.F90:37
pwadvection_mod::l_tend_pr_tot_w
logical l_tend_pr_tot_w
Definition: pwadvection.F90:37
pwadvection_mod::iqi
integer iqi
Definition: pwadvection.F90:41
pwadvection_mod::l_toplevel
logical l_toplevel
Definition: pwadvection.F90:21
pwadvection_mod::tend_pr_tot_th
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_th
Definition: pwadvection.F90:33
collections_mod
Collection data structures.
Definition: collections.F90:7
pwadvection_mod::pwadvection_get_descriptor
type(component_descriptor_type) function, public pwadvection_get_descriptor()
Provides a description of this component for the core to register.
Definition: pwadvection.F90:51
pwadvection_mod::tend_3d_u
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_u
Definition: pwadvection.F90:25
pwadvection_mod::l_tend_pr_tot_qg
logical l_tend_pr_tot_qg
Definition: pwadvection.F90:37
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
grids_mod::y_index
integer, parameter, public y_index
Definition: grids.F90:14
pwadvection_mod::l_tend_pr_tot_qi
logical l_tend_pr_tot_qi
Definition: pwadvection.F90:37
pwadvection_mod::determine_if_advection_here
logical function determine_if_advection_here(field)
Parses a field string (read in from the configuration file) and determines whether this algorithm sho...
Definition: pwadvection.F90:926
pwadvection_mod::save_precomponent_tendencies
subroutine save_precomponent_tendencies(current_state, cxn, cyn, txn, tyn)
Save the 3d tendencies coming into the component.
Definition: pwadvection.F90:625
pwadvection_mod::iqr
integer iqr
Definition: pwadvection.F90:41
pwadvection_mod::iql
integer iql
Definition: pwadvection.F90:41
pwadvection_mod::tend_3d_qv
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_qv
Definition: pwadvection.F90:25
pwadvection_mod::l_tend_3d_u
logical l_tend_3d_u
Definition: pwadvection.F90:29
pwadvection_mod::advect_flow
logical advect_flow
Definition: pwadvection.F90:19
pwadvection_mod::tend_3d_w
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_w
Definition: pwadvection.F90:25
pwadvection_mod::l_tend_3d_qv
logical l_tend_3d_qv
Definition: pwadvection.F90:29
pwadvection_mod::tend_3d_ql
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_ql
Definition: pwadvection.F90:25
pwadvection_mod::tend_pr_tot_ql
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_ql
Definition: pwadvection.F90:33
pwadvection_mod::tend_3d_qi
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_qi
Definition: pwadvection.F90:25
pwadvection_mod::l_tend_pr_tot_v
logical l_tend_pr_tot_v
Definition: pwadvection.F90:37
monc_component_mod
Interfaces and types that MONC components must specify.
Definition: monc_component.F90:6
pwadvection_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: pwadvection.F90:844
pwadvection_mod::advect_q_field
subroutine advect_q_field(current_state, current_x_index, current_y_index)
Advects the q fields in a column.
Definition: pwadvection.F90:333
pwadvection_mod::timestep_callback
subroutine timestep_callback(current_state)
Called per column of data, this will perform Piacsek-Williams advection on the applicable fields for ...
Definition: pwadvection.F90:264
pwadvection_mod::l_tend_3d_qs
logical l_tend_3d_qs
Definition: pwadvection.F90:29
pwadvection_mod::iqs
integer iqs
Definition: pwadvection.F90:41
optionsdatabase_mod::options_get_string
character(len=string_length) function, public options_get_string(options_database, key, index)
Retrieves a string value from the database that matches the provided key.
Definition: optionsdatabase.F90:280
pwadvection_mod::tend_3d_qs
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_qs
Definition: pwadvection.F90:25
pwadvection_mod::advect_q
logical advect_q
Definition: pwadvection.F90:19
pwadvection_mod::l_tend_pr_tot_ql
logical l_tend_pr_tot_ql
Definition: pwadvection.F90:37
pwadvection_mod::l_tend_pr_tot_qs
logical l_tend_pr_tot_qs
Definition: pwadvection.F90:37
pwadvection_mod::l_tend_pr_tot_tabs
logical l_tend_pr_tot_tabs
Definition: pwadvection.F90:37
pwadvection_mod::l_tend_3d_w
logical l_tend_3d_w
Definition: pwadvection.F90:29
pwadvection_mod::advect_th
logical advect_th
Definition: pwadvection.F90:19
pwadvection_mod::tend_pr_tot_qg
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_qg
Definition: pwadvection.F90:33
pwadvection_mod::diagnostic_generation_frequency
integer diagnostic_generation_frequency
Definition: pwadvection.F90:42
grids_mod::z_index
integer, parameter, public z_index
Grid index parameters.
Definition: grids.F90:14
pwadvection_mod::tend_pr_tot_qv
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_qv
Definition: pwadvection.F90:33
pwadvection_mod::tend_3d_qr
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_qr
Definition: pwadvection.F90:25
pwadvection_mod::tend_3d_th
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_th
Definition: pwadvection.F90:25
pwadvection_mod::l_tend_3d_qr
logical l_tend_3d_qr
Definition: pwadvection.F90:29
q_indices_mod::standard_q_names
type(standard_q_names_type), public standard_q_names
Definition: q_indices.F90:59
pwadvection_mod::l_tend_3d_th
logical l_tend_3d_th
Definition: pwadvection.F90:29
pwadvection_mod::compute_component_tendencies
subroutine compute_component_tendencies(current_state, cxn, cyn, txn, tyn)
Computation of component tendencies.
Definition: pwadvection.F90:673
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
monc_component_mod::component_field_information_type
Definition: monc_component.F90:31
pwadvection_mod::tend_pr_tot_w
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_w
Definition: pwadvection.F90:33
pwadvection_mod::l_tend_3d_qi
logical l_tend_3d_qi
Definition: pwadvection.F90:29
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
pwadvection_mod::l_tend_pr_tot_qr
logical l_tend_pr_tot_qr
Definition: pwadvection.F90:37
pwadvection_mod::finalisation_callback
subroutine finalisation_callback(current_state)
Definition: pwadvection.F90:230
datadefn_mod
Contains common definitions for the data and datatypes used by MONC.
Definition: datadefn.F90:2
pwadvection_mod::tend_pr_tot_qr
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_qr
Definition: pwadvection.F90:33
pwadvection_mod::l_tend_3d_tabs
logical l_tend_3d_tabs
Definition: pwadvection.F90:29
pwadvection_mod::advect_th_field
subroutine advect_th_field(current_state, current_x_index, current_y_index)
Advects the theta field in a column.
Definition: pwadvection.F90:403
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
pwadvection_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: pwadvection.F90:757
pwadvection_mod::l_tend_pr_tot_u
logical l_tend_pr_tot_u
Definition: pwadvection.F90:37
pwadvection_mod::iqv
integer iqv
Definition: pwadvection.F90:41
pwadvection_mod::tend_pr_tot_u
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_u
Definition: pwadvection.F90:33
pwadvection_mod::tend_3d_tabs
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_tabs
Definition: pwadvection.F90:25
pwadvection_mod::tend_pr_tot_qi
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_qi
Definition: pwadvection.F90:33
grids_mod
Functionality to support the different types of grid and abstraction between global grids and local o...
Definition: grids.F90:5
optionsdatabase_mod
Manages the options database. Contains administration functions and deduce runtime options from the c...
Definition: optionsdatabase.F90:7
pwadvection_mod::initialisation_callback
subroutine initialisation_callback(current_state)
Initialisation callback, will set up the configuration of this advection scheme.
Definition: pwadvection.F90:89
pwadvection_mod::l_tend_3d_qg
logical l_tend_3d_qg
Definition: pwadvection.F90:29
pwadvection_mod::tend_pr_tot_qs
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_qs
Definition: pwadvection.F90:33
pwadvection_mod::l_tend_3d_v
logical l_tend_3d_v
Definition: pwadvection.F90:29
pwadvection_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: pwadvection.F90:905
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
pwadvection_mod::tend_3d_v
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_v
Definition: pwadvection.F90:25
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
pwadvection_mod::tend_pr_tot_tabs
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_tabs
Definition: pwadvection.F90:33