MONC
coriolis.F90
Go to the documentation of this file.
1 
8  use state_mod, only : model_state_type
9  use grids_mod, only : z_index, y_index, x_index
10  implicit none
11 
12 #ifndef TEST_MODE
13  private
14 #endif
15 
17  real(kind=default_precision), dimension(:), allocatable :: geostrophic_wind_x, geostrophic_wind_y
18  real(kind=default_precision) :: fcoriol
19 
20  ! Local tendency diagnostic variables for this component
21  ! 3D tendency fields and logicals for their use
22  real(kind=default_precision), dimension(:,:,:), allocatable :: &
25  ! Local mean tendency profile fields and logicals for their use
26  real(kind=default_precision), dimension(:), allocatable :: &
29 
31 
33 
34  contains
35 
39  coriolis_get_descriptor%name="coriolis"
40  coriolis_get_descriptor%version=0.1
44 
47  allocate(coriolis_get_descriptor%published_fields(2+2))
48 
49  coriolis_get_descriptor%published_fields(1)= "tend_u_coriolis_3d_local"
50  coriolis_get_descriptor%published_fields(2)= "tend_v_coriolis_3d_local"
51 
52  coriolis_get_descriptor%published_fields(2+1)= "tend_u_coriolis_profile_total_local"
53  coriolis_get_descriptor%published_fields(2+2)= "tend_v_coriolis_profile_total_local"
54 
55  end function coriolis_get_descriptor
56 
57 
62  subroutine field_information_retrieval_callback(current_state, name, field_information)
63  type(model_state_type), target, intent(inout) :: current_state
64  character(len=*), intent(in) :: name
65  type(component_field_information_type), intent(out) :: field_information
66  integer :: strcomp
67 
68  ! Field information for 3d
69  strcomp=index(name, "coriolis_3d_local")
70  if (strcomp .ne. 0) then
71  field_information%field_type=component_array_field_type
72  field_information%number_dimensions=3
73  field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
74  field_information%dimension_sizes(2)=current_state%local_grid%size(y_index)
75  field_information%dimension_sizes(3)=current_state%local_grid%size(x_index)
76  field_information%data_type=component_double_data_type
77 
78  if (name .eq. "tend_u_coriolis_3d_local") then
79  field_information%enabled=l_tend_3d_u
80  else if (name .eq. "tend_v_coriolis_3d_local") then
81  field_information%enabled=l_tend_3d_v
82  else
83  field_information%enabled=.true.
84  end if
85 
86  end if !end 3d check
87 
88  ! Field information for profiles
89  strcomp=index(name, "coriolis_profile_total_local")
90  if (strcomp .ne. 0) then
91  field_information%field_type=component_array_field_type
92  field_information%number_dimensions=1
93  field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
94  field_information%data_type=component_double_data_type
95 
96  if (name .eq. "tend_u_coriolis_profile_total_local") then
97  field_information%enabled=l_tend_pr_tot_u
98  else if (name .eq. "tend_v_coriolis_profile_total_local") then
99  field_information%enabled=l_tend_pr_tot_v
100  else
101  field_information%enabled=.true.
102  end if
103 
104  end if !end profile check
105 
107 
108 
113  subroutine field_value_retrieval_callback(current_state, name, field_value)
114  type(model_state_type), target, intent(inout) :: current_state
115  character(len=*), intent(in) :: name
116  type(component_field_value_type), intent(out) :: field_value
117 
118  ! 3d Tendency Fields
119  if (name .eq. "tend_u_coriolis_3d_local" .and. allocated(tend_3d_u)) then
120  call set_published_field_value(field_value, real_3d_field=tend_3d_u)
121  else if (name .eq. "tend_v_coriolis_3d_local" .and. allocated(tend_3d_v)) then
122  call set_published_field_value(field_value, real_3d_field=tend_3d_v)
123 
124  ! Profile Tendency Fields
125  else if (name .eq. "tend_u_coriolis_profile_total_local" .and. allocated(tend_pr_tot_u)) then
126  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_u)
127  else if (name .eq. "tend_v_coriolis_profile_total_local" .and. allocated(tend_pr_tot_v)) then
128  call set_published_field_value(field_value, real_1d_field=tend_pr_tot_v)
129  end if
130 
131 
132  end subroutine field_value_retrieval_callback
133 
134 
137  subroutine initialisation_callback(current_state)
138  type(model_state_type), target, intent(inout) :: current_state
139 
140  integer :: k
141 
142  baroclinicity_use_geostrophic_shear=options_get_logical(current_state%options_database, "baroclinicity_use_geostrophic_shear")
143  fcoriol=options_get_real(current_state%options_database, "fcoriol")
144  current_state%geostrophic_wind_rate_of_change_in_x=options_get_real(current_state%options_database, &
145  "geostrophic_wind_rate_of_change_in_x")
146  current_state%geostrophic_wind_rate_of_change_in_y=options_get_real(current_state%options_database, &
147  "geostrophic_wind_rate_of_change_in_y")
148  current_state%surface_geostrophic_wind_x=options_get_real(current_state%options_database, "surface_geostrophic_wind_x")
149  current_state%surface_geostrophic_wind_y=options_get_real(current_state%options_database, "surface_geostrophic_wind_y")
150 
151  allocate(geostrophic_wind_x(current_state%local_grid%size(z_index)), &
152  geostrophic_wind_y(current_state%local_grid%size(z_index)))
153 
154  do k=1,current_state%local_grid%size(z_index)
155  geostrophic_wind_x(k)=current_state%surface_geostrophic_wind_x
156  geostrophic_wind_y(k)=current_state%surface_geostrophic_wind_y
158  geostrophic_wind_x(k)=geostrophic_wind_x(k)+current_state%geostrophic_wind_rate_of_change_in_x*&
159  current_state%global_grid%configuration%vertical%zn(k)
160  geostrophic_wind_y(k)=geostrophic_wind_y(k)+current_state%geostrophic_wind_rate_of_change_in_y*&
161  current_state%global_grid%configuration%vertical%zn(k)
162  end if
163  end do
164 
165  ! Tendency Logicals
166  l_tend_pr_tot_u = current_state%u%active
167  l_tend_pr_tot_v = current_state%v%active
168 
169  l_tend_3d_u = current_state%u%active .or. l_tend_pr_tot_u
170  l_tend_3d_v = current_state%v%active .or. l_tend_pr_tot_v
171 
172  ! Allocate 3d tendency fields upon availability
173  if (l_tend_3d_u) then
174  allocate( tend_3d_u(current_state%local_grid%size(z_index), &
175  current_state%local_grid%size(y_index), &
176  current_state%local_grid%size(x_index) ) )
177  endif
178  if (l_tend_3d_v) then
179  allocate( tend_3d_v(current_state%local_grid%size(z_index), &
180  current_state%local_grid%size(y_index), &
181  current_state%local_grid%size(x_index) ) )
182  endif
183 
184  ! Allocate profile tendency fields upon availability
185  if (l_tend_pr_tot_u) then
186  allocate( tend_pr_tot_u(current_state%local_grid%size(z_index)) )
187  endif
188  if (l_tend_pr_tot_v) then
189  allocate( tend_pr_tot_v(current_state%local_grid%size(z_index)) )
190  endif
191 
192  ! Save the sampling_frequency to force diagnostic calculation on select time steps
193  diagnostic_generation_frequency=options_get_integer(current_state%options_database, "sampling_frequency")
194 
195  end subroutine initialisation_callback
196 
197 
198  subroutine finalisation_callback(current_state)
199  type(model_state_type), target, intent(inout) :: current_state
200 
201  if (allocated(tend_3d_u)) deallocate(tend_3d_u)
202  if (allocated(tend_3d_v)) deallocate(tend_3d_v)
203 
204  if (allocated(tend_pr_tot_u)) deallocate(tend_pr_tot_u)
205  if (allocated(tend_pr_tot_v)) deallocate(tend_pr_tot_v)
206 
207  end subroutine finalisation_callback
208 
211  subroutine timestep_callback(current_state)
212  type(model_state_type), target, intent(inout) :: current_state
213 
214  integer :: local_y, locaL_x, k, target_x_index, target_y_index
215 
216  local_y=current_state%column_local_y
217  local_x=current_state%column_local_x
218  target_y_index=local_y-current_state%local_grid%halo_size(y_index)
219  target_x_index=local_x-current_state%local_grid%halo_size(x_index)
220 
221  ! Zero profile tendency totals on first instance in the sum
222  if (current_state%first_timestep_column) then
223  if (l_tend_pr_tot_u) then
224  tend_pr_tot_u(:)= 0.0_default_precision
225  endif
226  if (l_tend_pr_tot_v) then
227  tend_pr_tot_v(:)= 0.0_default_precision
228  endif
229  endif ! zero totals
230 
231  if (current_state%halo_column) then
232  if (.not. ((current_state%column_local_y == current_state%local_grid%halo_size(y_index) .and. &
233  current_state%column_local_x .le. current_state%local_grid%local_domain_end_index(x_index) .and. &
234  current_state%column_local_x .ge. current_state%local_grid%local_domain_start_index(x_index)-1) .or. &
235  (current_state%column_local_x == current_state%local_grid%halo_size(x_index) .and. &
236  current_state%column_local_y .ge. current_state%local_grid%local_domain_start_index(y_index) &
237  .and. current_state%column_local_y .le. current_state%local_grid%local_domain_end_index(y_index)) )) return
238  end if
239 
240  if (mod(current_state%timestep, diagnostic_generation_frequency) == 0 .and. .not. current_state%halo_column) then
241  call save_precomponent_tendencies(current_state, local_x, local_y, target_x_index, target_y_index)
242  end if
243 
244  do k=2,current_state%local_grid%size(z_index)
245 #if defined(U_ACTIVE) && defined(V_ACTIVE)
246  current_state%su%data(k, current_state%column_local_y, current_state%column_local_x)=&
247  current_state%su%data(k, current_state%column_local_y, current_state%column_local_x)+fcoriol*&
248  (0.25_default_precision*(current_state%v%data(k, current_state%column_local_y, current_state%column_local_x)+&
249  current_state%v%data(k, current_state%column_local_y, current_state%column_local_x+1)+&
250  current_state%v%data(k, current_state%column_local_y-1, current_state%column_local_x)+&
251  current_state%v%data(k, current_state%column_local_y-1, current_state%column_local_x+1))+current_state%vgal-&
253 
254  current_state%sv%data(k, current_state%column_local_y, current_state%column_local_x)=&
255  current_state%sv%data(k, current_state%column_local_y, current_state%column_local_x)-fcoriol*&
256  (0.25_default_precision*(current_state%u%data(k, current_state%column_local_y, current_state%column_local_x)+&
257  current_state%u%data(k, current_state%column_local_y, current_state%column_local_x-1)+&
258  current_state%u%data(k, current_state%column_local_y+1, current_state%column_local_x)+&
259  current_state%u%data(k, current_state%column_local_y+1, current_state%column_local_x-1))+current_state%ugal-&
261 #endif
262  end do
263 
264  if (mod(current_state%timestep, diagnostic_generation_frequency) == 0 .and. .not. current_state%halo_column) then
265  call compute_component_tendencies(current_state, local_x, local_y, target_x_index, target_y_index)
266  end if
267 
268  end subroutine timestep_callback
269 
270 
277  subroutine save_precomponent_tendencies(current_state, cxn, cyn, txn, tyn)
278  type(model_state_type), target, intent(in) :: current_state
279  integer, intent(in) :: cxn, cyn, txn, tyn
280 
281  ! Save 3d tendency fields upon request (of 3d or profiles) and availability
282  if (l_tend_3d_u) then
283  tend_3d_u(:,tyn,txn)=current_state%su%data(:,cyn,cxn)
284  endif
285  if (l_tend_3d_v) then
286  tend_3d_v(:,tyn,txn)=current_state%sv%data(:,cyn,cxn)
287  endif
288 
289  end subroutine save_precomponent_tendencies
290 
291 
298  subroutine compute_component_tendencies(current_state, cxn, cyn, txn, tyn)
299  type(model_state_type), target, intent(inout) :: current_state
300  integer, intent(in) :: cxn, cyn, txn, tyn
301 
302  ! Calculate change in tendency due to component
303  if (l_tend_3d_u) then
304  tend_3d_u(:,tyn,txn)=current_state%su%data(:,cyn,cxn) - tend_3d_u(:,tyn,txn)
305  endif
306  if (l_tend_3d_v) then
307  tend_3d_v(:,tyn,txn)=current_state%sv%data(:,cyn,cxn) - tend_3d_v(:,tyn,txn)
308  endif
309 
310  ! Add local tendency fields to the profile total
311  if (l_tend_pr_tot_u) then
312  tend_pr_tot_u(:)=tend_pr_tot_u(:) + tend_3d_u(:,tyn,txn)
313  endif
314  if (l_tend_pr_tot_v) then
315  tend_pr_tot_v(:)=tend_pr_tot_v(:) + tend_3d_v(:,tyn,txn)
316  endif
317 
318  end subroutine compute_component_tendencies
319 
324  subroutine set_published_field_value(field_value, real_1d_field, real_2d_field, real_3d_field)
325  type(component_field_value_type), intent(inout) :: field_value
326  real(kind=default_precision), dimension(:), optional :: real_1d_field
327  real(kind=default_precision), dimension(:,:), optional :: real_2d_field
328  real(kind=default_precision), dimension(:,:,:), optional :: real_3d_field
329 
330  if (present(real_1d_field)) then
331  allocate(field_value%real_1d_array(size(real_1d_field)), source=real_1d_field)
332  else if (present(real_2d_field)) then
333  allocate(field_value%real_2d_array(size(real_2d_field, 1), size(real_2d_field, 2)), source=real_2d_field)
334  else if (present(real_3d_field)) then
335  allocate(field_value%real_3d_array(size(real_3d_field, 1), size(real_3d_field, 2), size(real_3d_field, 3)), &
336  source=real_3d_field)
337  end if
338  end subroutine set_published_field_value
339 
340 
341 end module coriolis_mod
coriolis_mod::finalisation_callback
subroutine finalisation_callback(current_state)
Definition: coriolis.F90:199
coriolis_mod::geostrophic_wind_y
real(kind=default_precision), dimension(:), allocatable geostrophic_wind_y
Definition: coriolis.F90:17
coriolis_mod::timestep_callback
subroutine timestep_callback(current_state)
For each none halo cell this will calculate the coriolis terms for su and sv fields.
Definition: coriolis.F90:212
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
coriolis_mod::coriolis_get_descriptor
type(component_descriptor_type) function, public coriolis_get_descriptor()
Provides the descriptor back to the caller and is used in component registration.
Definition: coriolis.F90:39
grids_mod::y_index
integer, parameter, public y_index
Definition: grids.F90:14
coriolis_mod::l_tend_3d_u
logical l_tend_3d_u
Definition: coriolis.F90:24
monc_component_mod
Interfaces and types that MONC components must specify.
Definition: monc_component.F90:6
coriolis_mod::l_tend_pr_tot_u
logical l_tend_pr_tot_u
Definition: coriolis.F90:28
coriolis_mod::tend_3d_u
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_u
Definition: coriolis.F90:22
coriolis_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 components published fi...
Definition: coriolis.F90:63
coriolis_mod::l_tend_3d_v
logical l_tend_3d_v
Definition: coriolis.F90:24
coriolis_mod::tend_pr_tot_v
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_v
Definition: coriolis.F90:26
coriolis_mod::l_tend_pr_tot_v
logical l_tend_pr_tot_v
Definition: coriolis.F90:28
coriolis_mod
This calculates the coriolis and mean pressure gradient terms which impact su and sv fields.
Definition: coriolis.F90:2
grids_mod::z_index
integer, parameter, public z_index
Grid index parameters.
Definition: grids.F90:14
coriolis_mod::diagnostic_generation_frequency
integer diagnostic_generation_frequency
Definition: coriolis.F90:30
coriolis_mod::baroclinicity_use_geostrophic_shear
logical baroclinicity_use_geostrophic_shear
Definition: coriolis.F90:16
coriolis_mod::compute_component_tendencies
subroutine compute_component_tendencies(current_state, cxn, cyn, txn, tyn)
Computation of component tendencies.
Definition: coriolis.F90:299
coriolis_mod::save_precomponent_tendencies
subroutine save_precomponent_tendencies(current_state, cxn, cyn, txn, tyn)
Save the 3d tendencies coming into the component.
Definition: coriolis.F90:278
coriolis_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: coriolis.F90:114
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
optionsdatabase_mod::options_get_logical
logical function, public options_get_logical(options_database, key, index)
Retrieves a logical value from the database that matches the provided key.
Definition: optionsdatabase.F90:154
coriolis_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: coriolis.F90:325
coriolis_mod::tend_pr_tot_u
real(kind=default_precision), dimension(:), allocatable tend_pr_tot_u
Definition: coriolis.F90:26
coriolis_mod::fcoriol
real(kind=default_precision) fcoriol
Definition: coriolis.F90:18
datadefn_mod
Contains common definitions for the data and datatypes used by MONC.
Definition: datadefn.F90:2
coriolis_mod::tend_3d_v
real(kind=default_precision), dimension(:,:,:), allocatable tend_3d_v
Definition: coriolis.F90:22
coriolis_mod::initialisation_callback
subroutine initialisation_callback(current_state)
Initialisation call back which will read in the coriolis configuration and set up the geostrophic win...
Definition: coriolis.F90:138
coriolis_mod::geostrophic_wind_x
real(kind=default_precision), dimension(:), allocatable geostrophic_wind_x
Definition: coriolis.F90:17
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
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
optionsdatabase_mod::options_get_real
real(kind=default_precision) function, public options_get_real(options_database, key, index)
Retrieves a real value from the database that matches the provided key.
Definition: optionsdatabase.F90:91
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