64 type(model_state_type),
target,
intent(inout) :: current_state
65 character(len=*),
intent(in) :: name
66 type(component_field_information_type),
intent(out) :: field_information
70 field_information%field_type=component_array_field_type
71 field_information%data_type=component_double_data_type
72 field_information%number_dimensions=1
73 field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
74 field_information%enabled=.true.
77 strcomp=index(name,
"buoyancy_3d_local")
78 if (strcomp .ne. 0)
then
79 field_information%field_type=component_array_field_type
80 field_information%number_dimensions=3
81 field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
82 field_information%dimension_sizes(2)=current_state%local_grid%size(y_index)
83 field_information%dimension_sizes(3)=current_state%local_grid%size(x_index)
84 field_information%data_type=component_double_data_type
86 if (name .eq.
"tend_w_buoyancy_3d_local")
then
89 field_information%enabled=.true.
95 strcomp=index(name,
"buoyancy_profile_total_local")
96 if (strcomp .ne. 0)
then
97 field_information%field_type=component_array_field_type
98 field_information%number_dimensions=1
99 field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
100 field_information%data_type=component_double_data_type
102 if (name .eq.
"tend_w_buoyancy_profile_total_local")
then
105 field_information%enabled=.true.
118 type(model_state_type),
target,
intent(inout) :: current_state
119 character(len=*),
intent(in) :: name
120 type(component_field_value_type),
intent(out) :: field_value
122 if (name .eq.
"w_buoyancy")
then
126 else if (name .eq.
"tend_w_buoyancy_3d_local" .and.
allocated(
tend_3d_w))
then
130 else if (name .eq.
"tend_w_buoyancy_profile_total_local" .and.
allocated(
tend_pr_tot_w))
then
140 type(model_state_type),
target,
intent(inout) :: current_state
144 if (.not. current_state%passive_q .and. current_state%number_q_fields > 0)
then
145 if (.not.
allocated(current_state%cq))
then
146 allocate(current_state%cq(current_state%number_q_fields))
147 current_state%cq=0.0_default_precision
149 iqv = get_q_index(standard_q_names%VAPOUR,
'buoyancy')
150 current_state%cq(
iqv) = ratio_mol_wts-1.0
154 z_size=current_state%global_grid%size(z_index)
163 allocate(
tend_3d_w(current_state%local_grid%size(z_index), &
164 current_state%local_grid%size(y_index), &
165 current_state%local_grid%size(x_index) ) )
170 allocate(
tend_pr_tot_w(current_state%local_grid%size(z_index)) )
180 type(model_state_type),
target,
intent(inout) :: current_state
196 type(model_state_type),
target,
intent(inout) :: current_state
199 integer :: current_x_index, current_y_index, target_x_index, target_y_index
201 current_x_index=current_state%column_local_x
202 current_y_index=current_state%column_local_y
203 target_y_index=current_y_index-current_state%local_grid%halo_size(y_index)
204 target_x_index=current_x_index-current_state%local_grid%halo_size(x_index)
207 if (current_state%first_timestep_column)
then
219 if (.not. current_state%passive_th .and. current_state%th%active)
then
220 do k=2,current_state%local_grid%size(z_index)-1
221 w_buoyancy(k)=(0.5_default_precision*current_state%global_grid%configuration%vertical%buoy_co(k))*&
222 (current_state%th%data(k, current_state%column_local_y, current_state%column_local_x)&
223 +current_state%th%data(k+1, current_state%column_local_y, current_state%column_local_x))
224 current_state%sw%data(k, current_state%column_local_y, current_state%column_local_x)=&
225 current_state%sw%data(k, current_state%column_local_y, current_state%column_local_x)+
w_buoyancy(k)
228 if (.not. current_state%passive_q .and. current_state%number_q_fields .gt. 0)
then
229 if (current_state%use_anelastic_equations)
then
230 do n=1,current_state%number_q_fields
231 do k=2,current_state%local_grid%size(z_index)-1
232 current_state%sw%data(k, current_state%column_local_y, current_state%column_local_x)=&
233 current_state%sw%data(k, current_state%column_local_y, current_state%column_local_x)+&
234 (0.5_default_precision*current_state%global_grid%configuration%vertical%buoy_co(k))*&
235 current_state%cq(n)* (current_state%global_grid%configuration%vertical%thref(k)*&
236 current_state%q(n)%data(k, current_state%column_local_y, current_state%column_local_x)+&
237 current_state%global_grid%configuration%vertical%thref(k+1)*&
238 current_state%q(n)%data(k+1, current_state%column_local_y, current_state%column_local_x))
242 do n=1,current_state%number_q_fields
243 do k=2,current_state%local_grid%size(z_index)-1
244 current_state%sw%data(k, current_state%column_local_y, current_state%column_local_x)=&
245 current_state%sw%data(k, current_state%column_local_y, current_state%column_local_x)+&
247 (current_state%q(n)%data(k, current_state%column_local_y, current_state%column_local_x)+&
248 current_state%q(n)%data(k+1, current_state%column_local_y, current_state%column_local_x))
269 type(model_state_type),
target,
intent(in) :: current_state
270 integer,
intent(in) :: cxn, cyn, txn, tyn
274 tend_3d_w(:,tyn,txn)=current_state%sw%data(:,cyn,cxn)
287 type(model_state_type),
target,
intent(inout) :: current_state
288 integer,
intent(in) :: cxn, cyn, txn, tyn
308 type(component_field_value_type),
intent(inout) :: field_value
309 real(kind=default_precision),
dimension(:),
optional :: real_1d_field
310 real(kind=default_precision),
dimension(:,:),
optional :: real_2d_field
311 real(kind=default_precision),
dimension(:,:,:),
optional :: real_3d_field
313 if (
present(real_1d_field))
then
314 allocate(field_value%real_1d_array(
size(real_1d_field)), source=real_1d_field)
315 else if (
present(real_2d_field))
then
316 allocate(field_value%real_2d_array(
size(real_2d_field, 1),
size(real_2d_field, 2)), source=real_2d_field)
317 else if (
present(real_3d_field))
then
318 allocate(field_value%real_3d_array(
size(real_3d_field, 1),
size(real_3d_field, 2),
size(real_3d_field, 3)), &
319 source=real_3d_field)