MONC
viscosity.F90
Go to the documentation of this file.
1 
8  use grids_mod, only : z_index, x_index, y_index
14  implicit none
15 
16 #ifndef TEST_MODE
17  private
18 #endif
19 
20  real(kind=DEFAULT_PRECISION), dimension(:), allocatable :: u_viscosity, v_viscosity, w_viscosity
21 
23 
24 contains
25 
29  viscosity_get_descriptor%name="viscosity"
30  viscosity_get_descriptor%version=0.1
34 
37  allocate(viscosity_get_descriptor%published_fields(3))
38  viscosity_get_descriptor%published_fields(1)="u_viscosity"
39  viscosity_get_descriptor%published_fields(2)="v_viscosity"
40  viscosity_get_descriptor%published_fields(3)="w_viscosity"
41  end function viscosity_get_descriptor
42 
47  subroutine field_information_retrieval_callback(current_state, name, field_information)
48  type(model_state_type), target, intent(inout) :: current_state
49  character(len=*), intent(in) :: name
50  type(component_field_information_type), intent(out) :: field_information
51 
52  ! Field description is the same regardless of the specific field being retrieved
53  field_information%field_type=component_array_field_type
54  field_information%data_type=component_double_data_type
55  field_information%number_dimensions=1
56  field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
57  field_information%enabled=.true.
59 
64  subroutine field_value_retrieval_callback(current_state, name, field_value)
65  type(model_state_type), target, intent(inout) :: current_state
66  character(len=*), intent(in) :: name
67  type(component_field_value_type), intent(out) :: field_value
68 
69  if (name .eq. "u_viscosity") then
70  allocate(field_value%real_1d_array(size(u_viscosity)), source=u_viscosity)
71  else if (name .eq. "v_viscosity") then
72  allocate(field_value%real_1d_array(size(v_viscosity)), source=v_viscosity)
73  else if (name .eq. "w_viscosity") then
74  allocate(field_value%real_1d_array(size(w_viscosity)), source=w_viscosity)
75  end if
76  end subroutine field_value_retrieval_callback
77 
80  subroutine initialisation_callback(current_state)
81  type(model_state_type), target, intent(inout) :: current_state
82 
83  integer :: z_size, y_size, x_size
84 
85  z_size=current_state%local_grid%size(z_index) + current_state%local_grid%halo_size(z_index) * 2
86  y_size=current_state%local_grid%size(y_index) + current_state%local_grid%halo_size(y_index) * 2
87  x_size=current_state%local_grid%size(x_index) + current_state%local_grid%halo_size(x_index) * 2
88  allocate(current_state%vis_coefficient%data(z_size, y_size, x_size))
89 
90  z_size=current_state%global_grid%size(z_index)
91  allocate(u_viscosity(z_size), v_viscosity(z_size), w_viscosity(z_size))
92 
93  end subroutine initialisation_callback
94 
95  subroutine finalisation_callback(current_state)
96  type(model_state_type), target, intent(inout) :: current_state
97 
98  if (allocated(u_viscosity)) deallocate(u_viscosity)
99  if (allocated(v_viscosity)) deallocate(v_viscosity)
100  if (allocated(w_viscosity)) deallocate(w_viscosity)
101  end subroutine finalisation_callback
102 
105  subroutine timestep_callback(current_state)
106  type(model_state_type), target, intent(inout) :: current_state
107 
108  integer :: local_y, locaL_x, k
109  real(kind=DEFAULT_PRECISION), dimension(current_state%local_grid%size(Z_INDEX)) :: tau12, tau12_ym1, tau12m1, &
110  tau11, tau22, tau22_yp1, tau33, tau23_ym1, tau11p1, tau13, tau13m1, tau23
111 
112  if (.not. current_state%use_viscosity_and_diffusion .or. current_state%halo_column) return
113  if (current_state%viscosity_halo_swap_state%swap_in_progress) then
114  ! If there is a viscosity halo swap in progress then complete it
115  call complete_nonblocking_halo_swap(current_state, current_state%viscosity_halo_swap_state, &
117  end if
118 
119  local_y=current_state%column_local_y
120  local_x=current_state%column_local_x
121  if (current_state%field_stepping == forward_stepping) then
122  call calculate_tau(current_state, local_y, local_x, current_state%u, current_state%v, current_state%w, tau12, tau12_ym1, &
123  tau12m1, tau11, tau22, tau22_yp1, tau33, tau11p1, tau13, tau13m1, tau23, tau23_ym1)
124  else
125  call calculate_tau(current_state, local_y, local_x, current_state%zu, current_state%zv, current_state%zw, tau12, tau12_ym1, &
126  tau12m1, tau11, tau22, tau22_yp1, tau33, tau11p1, tau13, tau13m1, tau23, tau23_ym1)
127  end if
128  call calculate_viscous_sources(current_state, local_y, local_x, tau12, tau12_ym1, tau12m1, tau11, tau22, tau22_yp1, tau33, &
129  tau11p1, tau13, tau13m1, tau23, tau23_ym1)
130  end subroutine timestep_callback
131 
136  subroutine calculate_viscous_sources(current_state, local_y, local_x, tau12, tau12_ym1, tau12m1, &
137  tau11, tau22, tau22_yp1, tau33, tau11p1, tau13, tau13m1, tau23, tau23_ym1)
138  type(model_state_type), target, intent(inout) :: current_state
139  integer, intent(in) :: local_y, local_x
140  real(kind=DEFAULT_PRECISION), dimension(:), intent(in) :: tau12, tau12_ym1, tau12m1, tau11, tau22, tau22_yp1, tau33, &
141  tau11p1, tau13, tau13m1, tau23, tau23_ym1
142 
143  integer :: k
144 
145  do k=2, current_state%local_grid%size(z_index)
146 #ifdef U_ACTIVE
147  u_viscosity(k)=((tau11p1(k)-tau11(k))*current_state%global_grid%configuration%horizontal%cx+(tau12(k)-tau12_ym1(k))*&
148  current_state%global_grid%configuration%horizontal%cy+(tau13(k)-tau13(k-1))*&
149  current_state%global_grid%configuration%vertical%rdz(k))/current_state%global_grid%configuration%vertical%rhon(k)
150  current_state%su%data(k, local_y, local_x)=current_state%su%data(k, local_y, local_x)+u_viscosity(k)
151 #endif
152 #ifdef V_ACTIVE
153  v_viscosity(k)=((tau12(k)-tau12m1(k))*current_state%global_grid%configuration%horizontal%cx+(tau22_yp1(k)-tau22(k))*&
154  current_state%global_grid%configuration%horizontal%cy+(tau23(k)-tau23(k-1))*&
155  current_state%global_grid%configuration%vertical%rdz(k))/current_state%global_grid%configuration%vertical%rhon(k)
156  current_state%sv%data(k, local_y, local_x)=current_state%sv%data(k, local_y, local_x)+v_viscosity(k)
157 #endif
158  end do
159 #ifdef W_ACTIVE
160  do k=2, current_state%local_grid%size(z_index)-1
161  w_viscosity(k)=((tau13(k)-tau13m1(k))*current_state%global_grid%configuration%horizontal%cx+(tau23(k)-tau23_ym1(k))*&
162  current_state%global_grid%configuration%horizontal%cy+(tau33(k+1)-tau33(k))*&
163  current_state%global_grid%configuration%vertical%rdzn(k+1))/current_state%global_grid%configuration%vertical%rho(k)
164  current_state%sw%data(k, local_y, local_x)=current_state%sw%data(k, local_y, local_x)+w_viscosity(k)
165  end do
166 #endif
167  end subroutine calculate_viscous_sources
168 
173  subroutine calculate_tau(current_state, local_y, local_x, zu, zv, zw, tau12, tau12_ym1, tau12m1, tau11, &
174  tau22, tau22_yp1, tau33, tau11p1, tau13, tau13m1, tau23, tau23_ym1)
175  type(model_state_type), target, intent(inout) :: current_state
176  type(prognostic_field_type), intent(inout) :: zu, zv, zw
177  integer, intent(in) :: local_y, local_x
178  real(kind=DEFAULT_PRECISION), dimension(:), intent(out) :: tau12, tau12m1, tau11, tau22, tau22_yp1, tau33, &
179  tau11p1, tau13, tau13m1, tau23, tau23_ym1, tau12_ym1
180 
181  integer :: k
182  real(kind=DEFAULT_PRECISION) :: vistmp, vis12, vis12a, visonp2, visonp2a, vis13, vis23
183 
184  ! Do p levels and w-levels
185  do k=1, current_state%local_grid%size(z_index)
186  if (k .gt. 1) then
187  vistmp=current_state%vis_coefficient%data(k, local_y, local_x)+current_state%vis_coefficient%data(k, local_y+1, local_x)+&
188  current_state%vis_coefficient%data(k-1, local_y, local_x)+current_state%vis_coefficient%data(k-1, local_y+1, local_x)
189  vis12=0.125_default_precision*(vistmp + current_state%vis_coefficient%data(k, local_y, local_x+1)+&
190  current_state%vis_coefficient%data(k, local_y+1, local_x+1)+&
191  current_state%vis_coefficient%data(k-1, local_y, local_x+1)+&
192  current_state%vis_coefficient%data(k-1, local_y+1, local_x+1))
193  tau12(k)=0.0_default_precision
194 #ifdef U_ACTIVE
195  tau12(k)=(zu%data(k, local_y+1, local_x)-zu%data(k, local_y, local_x))*&
196  current_state%global_grid%configuration%horizontal%cy
197 #endif
198 #ifdef V_ACTIVE
199  tau12(k)=tau12(k)+(zv%data(k, local_y, local_x+1)-zv%data(k, local_y, local_x))*&
200  current_state%global_grid%configuration%horizontal%cx
201 #endif
202  tau12(k)=tau12(k)*current_state%global_grid%configuration%vertical%rhon(k)*vis12
203 
204  vis12a=0.125_default_precision*(vistmp +current_state%vis_coefficient%data(k, local_y, local_x-1)+&
205  current_state%vis_coefficient%data(k, local_y+1, local_x-1)+&
206  current_state%vis_coefficient%data(k-1, local_y, local_x-1)+&
207  current_state%vis_coefficient%data(k-1, local_y+1, local_x-1))
208  tau12m1(k)=0.0_default_precision
209 #ifdef U_ACTIVE
210  tau12m1(k)=(zu%data(k, local_y+1, local_x-1)-zu%data(k, local_y, local_x-1))*&
211  current_state%global_grid%configuration%horizontal%cy
212 #endif
213 #ifdef V_ACTIVE
214  tau12m1(k)=tau12m1(k)+(zv%data(k, local_y, local_x)-zv%data(k, local_y, local_x-1))*&
215  current_state%global_grid%configuration%horizontal%cx
216 #endif
217  tau12m1(k)=tau12m1(k)*current_state%global_grid%configuration%vertical%rhon(k)*vis12a
218 
219  vistmp=current_state%vis_coefficient%data(k, local_y-1, local_x)+current_state%vis_coefficient%data(k, local_y, local_x)+&
220  current_state%vis_coefficient%data(k-1, local_y-1, local_x)+current_state%vis_coefficient%data(k-1, local_y, local_x)
221  vis12=0.125_default_precision*(vistmp + current_state%vis_coefficient%data(k, local_y-1, local_x+1)+&
222  current_state%vis_coefficient%data(k, local_y, local_x+1)+&
223  current_state%vis_coefficient%data(k-1, local_y-1, local_x+1)+&
224  current_state%vis_coefficient%data(k-1, local_y, local_x+1))
225  tau12_ym1(k)=0.0_default_precision
226 #ifdef U_ACTIVE
227  tau12_ym1(k)=(zu%data(k, local_y, local_x)-zu%data(k, local_y-1, local_x))*&
228  current_state%global_grid%configuration%horizontal%cy
229 #endif
230 #ifdef V_ACTIVE
231  tau12_ym1(k)=tau12_ym1(k)+(zv%data(k, local_y-1, local_x+1)-zv%data(k, local_y-1, local_x))*&
232  current_state%global_grid%configuration%horizontal%cx
233 #endif
234  tau12_ym1(k)=tau12_ym1(k)*current_state%global_grid%configuration%vertical%rhon(k)*vis12
235 
236  visonp2=current_state%global_grid%configuration%vertical%rhon(k)*&
237  (current_state%vis_coefficient%data(k, local_y, local_x)+current_state%vis_coefficient%data(k-1, local_y, local_x))
238 #ifdef U_ACTIVE
239  tau11(k)=visonp2*(zu%data(k, local_y, local_x)-zu%data(k, local_y, local_x-1))*&
240  current_state%global_grid%configuration%horizontal%cx
241 #else
242  tau11(k)=0.0_default_precision
243 #endif
244 #ifdef V_ACTIVE
245  tau22(k)=visonp2*(zv%data(k, local_y, local_x)-zv%data(k, local_y-1, local_x))*&
246  current_state%global_grid%configuration%horizontal%cy
247 #else
248  tau22(k)=0.0_default_precision
249 #endif
250 #ifdef W_ACTIVE
251  tau33(k)=visonp2*(zw%data(k, local_y, local_x)-zw%data(k-1, local_y, local_x))*&
252  current_state%global_grid%configuration%vertical%rdz(k)
253 #else
254  tau33(k)=0.0_default_precision
255 #endif
256 
257 #ifdef V_ACTIVE
258  visonp2=current_state%global_grid%configuration%vertical%rhon(k)*&
259  (current_state%vis_coefficient%data(k, local_y+1, local_x)+&
260  current_state%vis_coefficient%data(k-1, local_y+1, local_x))
261  tau22_yp1(k)=visonp2*(zv%data(k, local_y+1, local_x)-zv%data(k, local_y, local_x))*&
262  current_state%global_grid%configuration%horizontal%cy
263 #endif
264 
265 #ifdef U_ACTIVE
266  visonp2a=current_state%global_grid%configuration%vertical%rhon(k)*&
267  (current_state%vis_coefficient%data(k, local_y, local_x+1)+&
268  current_state%vis_coefficient%data(k-1, local_y, local_x+1))
269  tau11p1(k)=visonp2a*(current_state%zu%data(k, local_y, local_x+1)-&
270  zu%data(k, local_y, local_x))*current_state%global_grid%configuration%horizontal%cx
271 #endif
272  else
273  tau12(k)=0.0_default_precision
274  tau12_ym1(k)=0.0_default_precision
275  tau12m1(k)=0.0_default_precision
276  tau11(k)=0.0_default_precision
277  tau22(k)=0.0_default_precision
278  tau22_yp1(k)=0.0_default_precision
279  tau33(k)=0.0_default_precision
280  tau11p1(k)=0.0_default_precision
281  end if
282  if (k .lt. current_state%local_grid%size(z_index)) then
283  vis13=0.5_default_precision*(current_state%vis_coefficient%data(k, local_y, local_x)+&
284  current_state%vis_coefficient%data(k, local_y, local_x+1))
285  tau13(k)=0.0_default_precision
286 #ifdef U_ACTIVE
287  tau13(k)=(zu%data(k+1, local_y, local_x)-zu%data(k, local_y, local_x))*&
288  current_state%global_grid%configuration%vertical%rdzn(k+1)
289 #endif
290 #ifdef W_ACTIVE
291  tau13(k)=tau13(k)+(zw%data(k, local_y, local_x+1)-zw%data(k, local_y, local_x))*&
292  current_state%global_grid%configuration%horizontal%cx
293 #endif
294  tau13(k)=tau13(k)*current_state%global_grid%configuration%vertical%rho(k)*vis13
295 
296  vis13=0.5_default_precision*(current_state%vis_coefficient%data(k, local_y, local_x-1)+&
297  current_state%vis_coefficient%data(k, local_y, local_x))
298  tau13m1(k)=0.0_default_precision
299 #ifdef U_ACTIVE
300  tau13m1(k)=(zu%data(k+1, local_y, local_x-1)-zu%data(k, local_y, local_x-1))*&
301  current_state%global_grid%configuration%vertical%rdzn(k+1)
302 #endif
303 #ifdef W_ACTIVE
304  tau13m1(k)=tau13m1(k)+(zw%data(k, local_y, local_x)-zw%data(k, local_y, local_x-1))*&
305  current_state%global_grid%configuration%horizontal%cx
306 #endif
307  tau13m1(k)=tau13m1(k)*current_state%global_grid%configuration%vertical%rho(k)*vis13
308 
309  vis23=0.5_default_precision*(current_state%vis_coefficient%data(k, local_y, local_x)+&
310  current_state%vis_coefficient%data(k, local_y+1, local_x))
311  tau23(k)=0.0_default_precision
312 #ifdef W_ACTIVE
313  tau23(k)=(zw%data(k, local_y+1, local_x)-zw%data(k, local_y, local_x))*&
314  current_state%global_grid%configuration%horizontal%cy
315 #endif
316 #ifdef V_ACTIVE
317  tau23(k)=tau23(k)+(zv%data(k+1, local_y, local_x)-zv%data(k, local_y, local_x))*&
318  current_state%global_grid%configuration%vertical%rdzn(k+1)
319 #endif
320  tau23(k)=tau23(k)*current_state%global_grid%configuration%vertical%rho(k)*vis23
321 
322  vis23=0.5_default_precision*(current_state%vis_coefficient%data(k, local_y-1, local_x)+&
323  current_state%vis_coefficient%data(k, local_y, local_x))
324  tau23_ym1(k)=0.0_default_precision
325 #ifdef W_ACTIVE
326  tau23_ym1(k)=(zw%data(k, local_y, local_x)-zw%data(k, local_y-1, local_x))*&
327  current_state%global_grid%configuration%horizontal%cy
328 #endif
329 #ifdef V_ACTIVE
330  tau23_ym1(k)=tau23_ym1(k)+(zv%data(k+1, local_y-1, local_x)-zv%data(k, local_y-1, local_x))*&
331  current_state%global_grid%configuration%vertical%rdzn(k+1)
332 #endif
333  tau23_ym1(k)=tau23_ym1(k)*current_state%global_grid%configuration%vertical%rho(k)*vis23
334  else
335  tau13(k)=0.0_default_precision
336  tau13m1(k)=0.0_default_precision
337  tau23(k)=0.0_default_precision
338  tau23_ym1(k)=0.0_default_precision
339  end if
340  end do
341  end subroutine calculate_tau
342 
346  subroutine perform_local_data_copy_for_vis(current_state, halo_depth, involve_corners, source_data)
347  type(model_state_type), intent(inout) :: current_state
348  integer, intent(in) :: halo_depth
349  logical, intent(in) :: involve_corners
350  type(field_data_wrapper_type), dimension(:), intent(in), optional :: source_data
351 
352  call perform_local_data_copy_for_field(current_state%vis_coefficient%data, current_state%local_grid, &
353  current_state%parallel%my_rank, halo_depth, involve_corners)
354  end subroutine perform_local_data_copy_for_vis
355 
364  subroutine copy_halo_buffer_to_vis(current_state, neighbour_description, dim, target_index, &
365  neighbour_location, current_page, source_data)
366  type(model_state_type), intent(inout) :: current_state
367  integer, intent(in) :: dim, target_index, neighbour_location
368  integer, intent(inout) :: current_page(:)
369  type(neighbour_description_type), intent(inout) :: neighbour_description
370  type(field_data_wrapper_type), dimension(:), intent(in), optional :: source_data
371 
372  call copy_buffer_to_field(current_state%local_grid, neighbour_description%recv_halo_buffer, &
373  current_state%vis_coefficient%data, dim, target_index, current_page(neighbour_location))
374 
375  current_page(neighbour_location)=current_page(neighbour_location)+1
376  end subroutine copy_halo_buffer_to_vis
377 
387  subroutine copy_halo_buffer_to_vis_corners(current_state, neighbour_description, corner_loc, x_target_index, &
388  y_target_index, neighbour_location, current_page, source_data)
389  type(model_state_type), intent(inout) :: current_state
390  integer, intent(in) :: corner_loc, x_target_index, y_target_index, neighbour_location
391  integer, intent(inout) :: current_page(:)
392  type(neighbour_description_type), intent(inout) :: neighbour_description
393  type(field_data_wrapper_type), dimension(:), intent(in), optional :: source_data
394 
395  call copy_buffer_to_corner(current_state%local_grid, neighbour_description%recv_corner_buffer, &
396  current_state%vis_coefficient%data, corner_loc, x_target_index, y_target_index, current_page(neighbour_location))
397 
398  current_page(neighbour_location)=current_page(neighbour_location)+1
399  end subroutine copy_halo_buffer_to_vis_corners
400 end module viscosity_mod
subroutine finalisation_callback(current_state)
Definition: viscosity.F90:96
Wrapper type for the value returned for a published field from a component.
subroutine calculate_viscous_sources(current_state, local_y, local_x, tau12, tau12_ym1, tau12m1, tau11, tau22, tau22_yp1, tau33, tau11p1, tau13, tau13m1, tau23, tau23_ym1)
Calculates the viscous sources based upon TAU for the U,V and W fields.
Definition: viscosity.F90:138
integer, parameter, public forward_stepping
Definition: state.F90:15
subroutine, public perform_local_data_copy_for_field(field_data, local_grid, my_rank, halo_depth, involve_corners)
Will perform a a local copy for the halo data of a field.
integer, parameter, public log_error
Only log ERROR messages.
Definition: logging.F90:11
Contains the types used for communication, holding the state of communications and supporting activit...
Contains prognostic field definitions and functions.
Definition: prognostics.F90:2
subroutine field_information_retrieval_callback(current_state, name, field_information)
Field information retrieval callback, this returns information for a specific components published fi...
Definition: viscosity.F90:48
A prognostic field which is assumed to be 3D.
Definition: prognostics.F90:13
Logging utility.
Definition: logging.F90:2
integer, parameter, public default_precision
MPI communication type which we use for the prognostic and calculation data.
Definition: datadefn.F90:17
integer, parameter, public z_index
Grid index parameters.
Definition: grids.F90:14
subroutine, public complete_nonblocking_halo_swap(current_state, halo_swap_state, perform_local_data_copy, copy_from_halo_buffer, copy_from_halo_buffer_to_corner, source_data)
This completes a nonblocking halo swap and it is only during this call that the data fields themselve...
type(component_descriptor_type) function, public viscosity_get_descriptor()
Provides the descriptor back to the caller and is used in component registration. ...
Definition: viscosity.F90:29
Contains common definitions for the data and datatypes used by MONC.
Definition: datadefn.F90:2
The ModelState which represents the current state of a run.
Definition: state.F90:39
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
subroutine, public log_master_log(level, message)
Will log just from the master process.
Definition: logging.F90:47
subroutine, public copy_buffer_to_corner(local_grid, halo_buffer, field_data, corner_loc, x_target_index, y_target_index, halo_page)
Copies the received buffer for a specific field to the corresponding corner of that field...
subroutine, public copy_buffer_to_field(local_grid, halo_buffer, field_data, dim, target_index, halo_page)
Copies the received buffer for a specific field to the corresponding halo data of that prognostic fie...
integer, parameter, public component_array_field_type
real(kind=default_precision), dimension(:), allocatable v_viscosity
Definition: viscosity.F90:20
Maintains the state of a halo swap and contains buffers, neighbours etc.
Interfaces and types that MONC components must specify.
real(kind=default_precision), dimension(:), allocatable w_viscosity
Definition: viscosity.F90:20
real(kind=default_precision), dimension(:), allocatable u_viscosity
Definition: viscosity.F90:20
Provides the mechanism for halo swapping. This module contains the functionality required to determin...
integer, parameter, public log_warn
Log WARNING and ERROR messages.
Definition: logging.F90:12
Describes the neighbours of a process in a specific dimension and contains the communication buffers ...
subroutine calculate_tau(current_state, local_y, local_x, zu, zv, zw, tau12, tau12_ym1, tau12m1, tau11, tau22, tau22_yp1, tau33, tau11p1, tau13, tau13m1, tau23, tau23_ym1)
Calculated TAU which is used in computing the viscous source terms.
Definition: viscosity.F90:175
subroutine perform_local_data_copy_for_vis(current_state, halo_depth, involve_corners, source_data)
Does local data copying for viscosity coefficient variable halo swap.
Definition: viscosity.F90:347
subroutine copy_halo_buffer_to_vis(current_state, neighbour_description, dim, target_index, neighbour_location, current_page, source_data)
Copies the halo buffer to halo location for the viscosity coefficient field.
Definition: viscosity.F90:366
Functionality to support the different types of grid and abstraction between global grids and local o...
Definition: grids.F90:5
subroutine copy_halo_buffer_to_vis_corners(current_state, neighbour_description, corner_loc, x_target_index, y_target_index, neighbour_location, current_page, source_data)
Copies the corner halo buffer to the viscosity coefficient field corners.
Definition: viscosity.F90:389
Computes the viscosity dynamics for the U,V,W source terms.
Definition: viscosity.F90:2
subroutine initialisation_callback(current_state)
Sets up the stencil_mod (used in interpolation) and allocates data for the flux fields.
Definition: viscosity.F90:81
subroutine timestep_callback(current_state)
At each timestep will compute the viscosity U,V,W source terms.
Definition: viscosity.F90:106
The model state which represents the current state of a run.
Definition: state.F90:2
integer, parameter, public y_index
Definition: grids.F90:14
integer, parameter, public x_index
Definition: grids.F90:14
MONC component registry.
Definition: registry.F90:5
integer, parameter, public component_double_data_type
subroutine field_value_retrieval_callback(current_state, name, field_value)
Field value retrieval callback, this returns the value of a specific published field.
Definition: viscosity.F90:65