11 use netcdf,
only : nf90_noerr, nf90_global, nf90_nowrite, nf90_inquire_attribute, nf90_open, nf90_strerror, &
12 nf90_inq_dimid, nf90_inquire_dimension, nf90_inq_varid, nf90_get_var, nf90_inquire, nf90_close, nf90_get_att
24 character(len=*),
parameter :: &
25 time_key =
"time", & !< NetCDF data time key
63 type(model_state_type),
intent(inout),
target :: current_state
68 if (.not. current_state%passive_q .and. current_state%number_q_fields > 0)
then
69 if (.not.
allocated(current_state%cq))
then
70 allocate(current_state%cq(current_state%number_q_fields))
71 current_state%cq=0.0_default_precision
73 iqv = get_q_index(standard_q_names%VAPOUR,
'setfluxlook')
74 current_state%cq(
iqv) = ratio_mol_wts-1.0
77 if (current_state%use_surface_boundary_conditions)
then
78 if (current_state%type_of_surface_boundary_conditions == prescribed_surface_fluxes)
then
79 if (current_state%use_time_varying_surface_values)
then
83 /(current_state%global_grid%configuration%vertical%rho(1)*cp)
85 /(current_state%global_grid%configuration%vertical%rho(1)*rlvap)
88 current_state%fbuoy=0.
89 if(.not. current_state%passive_th) current_state%fbuoy=&
90 current_state%global_grid%configuration%vertical%buoy_co(1)*current_state%surface_temperature_flux
91 if(.not. current_state%passive_q .and. current_state%number_q_fields > 0)
then
92 current_state%fbuoy=current_state%fbuoy+current_state%cq(
iqv)*current_state%surface_vapour_flux*g
95 current_state%theta_surf=0.0_default_precision
96 current_state%surface_vapour_mixing_ratio=0.0_default_precision
98 else if (current_state%type_of_surface_boundary_conditions == prescribed_surface_values)
then
101 if (current_state%use_time_varying_surface_values)
then
107 if (current_state%saturated_surface)
then
108 current_state%surface_vapour_mixing_ratio = qsaturation(
surface_temperatures(1),current_state%surface_pressure*0.01)
110 current_state%surface_vapour_mixing_ratio = &
111 options_get_real(current_state%options_database,
"surface_vapour_mixing_ratio")
118 (current_state%surface_reference_pressure/current_state%surface_pressure)**r_over_cp
119 current_state%theta_virtual_surf = current_state%theta_surf
120 if (.not. current_state%passive_q .and. current_state%number_q_fields > 0)
then
121 current_state%theta_virtual_surf = current_state%theta_surf + &
122 current_state%global_grid%configuration%vertical%thref(2)* &
123 current_state%cq(
iqv)*current_state%surface_vapour_mixing_ratio
127 current_state%cmbc=betam*current_state%global_grid%configuration%vertical%zn(2)*g*&
128 von_karman_constant/current_state%theta_virtual_surf
129 current_state%rcmbc=1.0_default_precision/current_state%cmbc
130 current_state%ellmocon=current_state%theta_virtual_surf/(g*von_karman_constant)
138 type(model_state_type),
intent(inout),
target :: current_state
140 if (current_state%use_surface_boundary_conditions)
then
141 if (current_state%use_time_varying_surface_values)
then
150 type(model_state_type),
intent(inout),
target :: current_state
152 integer I, ik, iters(current_state%lookup_table_entries)
153 real(kind=default_precision) :: smth, &
156 if (current_state%fbuoy .le. 0.0_default_precision)
return
157 current_state%velmax=100.0_default_precision
158 current_state%velmin=0.1_default_precision
159 current_state%aloginv=1.0_default_precision/log(current_state%velmax/current_state%velmin)
160 smth=0.1_default_precision
161 do ik=1, current_state%lookup_table_entries
162 current_state%lookup_table_velocity(ik)=current_state%velmin*(current_state%velmax/current_state%velmin)**&
163 (real(ik-1)/real(current_state%lookup_table_entries-1))
164 current_state%lookup_table_ustr(ik)=current_state%lookup_table_velocity(ik)*&
165 current_state%global_grid%configuration%vertical%vk_on_zlogm
166 if (current_state%fbuoy .gt. 0.0_default_precision)
then
169 iters(ik)=iters(ik)+1
170 ob=-current_state%lookup_table_ustr(ik)**3/(von_karman_constant*current_state%fbuoy)
171 x1=sqrt(sqrt(1.-gammam*(current_state%global_grid%configuration%vertical%zn(2)+z0)/ob))
172 x0=sqrt(sqrt(1.-gammam*z0/ob))
173 current_state%lookup_table_ustr(ik)=(1.-smth)*current_state%lookup_table_ustr(ik) + smth*&
174 (von_karman_constant*current_state%lookup_table_velocity(ik)) / &
175 (current_state%global_grid%configuration%vertical%zlogm-(2.0_default_precision*&
176 log((x1+1.0_default_precision)/(x0+1.0_default_precision)) + log((x1*x1+1.0_default_precision)/&
177 (x0*x0+1.0_default_precision)) + 2.0_default_precision*atan(x0) -2.0_default_precision*atan(x1)))
181 current_state%cneut=current_state%lookup_table_ustr(current_state%lookup_table_entries)/&
182 current_state%lookup_table_velocity(current_state%lookup_table_entries)
183 current_state%cfc=current_state%lookup_table_ustr(1)*current_state%lookup_table_velocity(1)**convective_limit
187 type(model_state_type),
intent(inout),
target :: current_state
189 real(kind=default_precision):: crelax,&
193 dvelustr(current_state%lookup_table_entries), ob, x1, x0
196 if (current_state%fbuoynew .le. 0.0_default_precision)
then
197 current_state%fbuoy=current_state%fbuoynew
206 if (nit .gt. 10 .or. (current_state%fbuoynew .gt. 0.0_default_precision .and. &
207 current_state%fbuoy .le. 0.0_default_precision))
then
208 current_state%fbuoy=current_state%fbuoynew
212 crelax=1./sqrt(real(nit))
213 dfb=(current_state%fbuoynew-current_state%fbuoy)/real(nit)
215 current_state%fbuoy=current_state%fbuoy+dfb
216 do ik=2, current_state%lookup_table_entries-1
217 dvelustr(ik)=log(current_state%lookup_table_velocity(ik+1)/current_state%lookup_table_velocity(ik-1))&
218 /log(current_state%lookup_table_ustr(ik+1)/current_state%lookup_table_ustr(ik-1))
220 dvelustr(1)=log(current_state%lookup_table_velocity(2)/current_state%lookup_table_velocity(1))/&
221 log(current_state%lookup_table_ustr(2)/current_state%lookup_table_ustr(1))
222 dvelustr(current_state%lookup_table_entries)=log(current_state%lookup_table_velocity(current_state%lookup_table_entries)/&
223 current_state%lookup_table_velocity(current_state%lookup_table_entries-1))/&
224 log(current_state%lookup_table_ustr(current_state%lookup_table_entries)/&
225 current_state%lookup_table_ustr(current_state%lookup_table_entries-1))
226 do ik=1, current_state%lookup_table_entries
228 ob=-current_state%lookup_table_ustr(ik)**3/(von_karman_constant*current_state%fbuoy)
229 x1=sqrt(sqrt(1.-gammam*(current_state%global_grid%configuration%vertical%zn(2)+z0)/ob))
230 x0=sqrt(sqrt(1.-gammam*z0/ob))
231 velnew=(current_state%lookup_table_ustr(ik)/von_karman_constant)*(current_state%global_grid%configuration%vertical%zlogm-&
232 (2.0_default_precision*log((x1+1.0_default_precision)/(x0+1.0_default_precision)) + &
233 log((x1*x1+1.0_default_precision)/(x0*x0+1.0_default_precision)) + 2.0_default_precision*atan(x0) &
234 -2.0_default_precision*atan(x1)))
236 current_state%lookup_table_ustr(ik)=current_state%lookup_table_ustr(ik)/&
237 (velnew/current_state%lookup_table_velocity(ik))**(crelax/dvelustr(ik))
240 current_state%cneut=current_state%lookup_table_ustr(current_state%lookup_table_entries)/&
241 current_state%lookup_table_velocity(current_state%lookup_table_entries)
242 current_state%cfc=current_state%lookup_table_ustr(1)*current_state%lookup_table_velocity(1)**convective_limit
243 current_state%fbuoy=current_state%fbuoynew
247 type(model_state_type),
intent(inout),
target :: current_state
249 real(kind=default_precision) :: surface_temp
251 if (current_state%type_of_surface_boundary_conditions == prescribed_surface_fluxes)
then
256 current_state%time, current_state%surface_temperature_flux, &
257 extrapolate=
'constant')
261 current_state%time, current_state%surface_vapour_flux, &
262 extrapolate=
'constant')
265 current_state%fbuoynew=0.0_default_precision
266 if (.not. current_state%passive_th) current_state%fbuoynew=&
267 current_state%global_grid%configuration%vertical%buoy_co(1)*current_state%surface_temperature_flux
268 if (.not. current_state%passive_q)
then
269 current_state%fbuoynew=current_state%fbuoynew+current_state%cq(
iqv)*current_state%surface_vapour_flux*g
272 else if (current_state%type_of_surface_boundary_conditions == prescribed_surface_values)
then
277 current_state%time, surface_temp, &
278 extrapolate=
'constant')
282 surface_temp = surface_temp + 273.15_default_precision
286 if (current_state%saturated_surface)
then
287 current_state%surface_vapour_mixing_ratio = qsaturation(surface_temp,current_state%surface_pressure*0.01)
291 current_state%time, current_state%surface_vapour_mixing_ratio, &
292 extrapolate=
'constant')
296 current_state%theta_surf = surface_temp*&
297 (current_state%surface_reference_pressure/current_state%surface_pressure)**r_over_cp
298 current_state%theta_virtual_surf = current_state%theta_surf
299 if (.not. current_state%passive_q)
then
300 current_state%theta_virtual_surf = current_state%theta_surf + &
301 current_state%global_grid%configuration%vertical%thref(2)* &
302 current_state%cq(
iqv)*current_state%surface_vapour_mixing_ratio
306 current_state%cmbc=betam*current_state%global_grid%configuration%vertical%zn(2)*g*&
307 von_karman_constant/current_state%theta_virtual_surf
308 current_state%rcmbc=1.0_default_precision/current_state%cmbc
309 current_state%ellmocon=current_state%theta_virtual_surf/(g*von_karman_constant)
316 type(model_state_type),
intent(inout),
target :: current_state
320 integer :: ncid, time_dim
321 integer :: number_input_humidities
324 current_state%use_surface_boundary_conditions= &
325 options_get_logical(current_state%options_database,
"use_surface_boundary_conditions")
327 if (current_state%use_surface_boundary_conditions)
then
328 current_state%type_of_surface_boundary_conditions=options_get_integer(current_state%options_database, &
329 "type_of_surface_boundary_conditions")
330 current_state%use_time_varying_surface_values=options_get_logical(current_state%options_database, &
331 "use_time_varying_surface_values")
333 current_state%saturated_surface = .true.
335 input_file=options_get_string(current_state%options_database,
"surface_conditions_file")
338 if (current_state%use_time_varying_surface_values)
then
343 if (current_state%type_of_surface_boundary_conditions == prescribed_surface_fluxes)
then
351 number_input_humidities=0
352 else if (current_state%type_of_surface_boundary_conditions == prescribed_surface_values)
then
358 call options_get_real_array(current_state%options_database,
"surface_temperatures",
surface_temperatures)
359 units_surface_temp=options_get_string(current_state%options_database,
"units_surface_temp")
361 call options_get_real_array(current_state%options_database,
"surface_humidities",
surface_humidities)
362 number_input_humidities=options_get_array_size(current_state%options_database,
"surface_humidities")
366 if (log_get_logging_level() .ge. log_debug)
then
367 call log_master_log(log_debug,
"Reading in surface boundary conditions from:"//trim(
input_file) )
374 number_input_humidities=0
380 if (number_input_humidities>0)current_state%saturated_surface=.false.
386 allocate(current_state%lookup_table_velocity(current_state%lookup_table_entries), &
387 current_state%lookup_table_ustr(current_state%lookup_table_entries))
397 integer,
intent(in) :: ncid
398 integer,
intent(out) :: time_dim
399 integer :: time_dimid
403 call check_status(nf90_inquire_dimension(ncid, time_dimid, len=time_dim))
412 subroutine read_variables(filename, ncid, time_dim, time, surface_temperatures, surface_humidities, &
413 surface_latent_heat_flux, surface_sensible_heat_flux )
414 character(*),
intent(in) :: filename
415 integer,
intent(in) :: ncid, time_dim
416 real(kind=default_precision),
dimension(:),
allocatable,
intent(inout) :: time
418 real(kind=default_precision),
dimension(:),
allocatable,
intent(inout) ::
surface_humidities
422 integer :: status, variable_id
428 status=nf90_inq_varid(ncid,
time_key, variable_id)
429 if (status==nf90_noerr)
then
430 allocate(time(time_dim))
433 call log_log(log_error,
"No recognized time variable found in"//trim(filename))
437 if (status==nf90_noerr)
then
443 if (status==nf90_noerr)
then
449 if (status==nf90_noerr)
then
455 if (status==nf90_noerr)
then
468 integer,
intent(in) :: ncid
469 character(len=*),
intent(in) :: key
470 real(kind=default_precision),
dimension(:),
intent(inout),
optional :: data1d
471 real(kind=default_precision),
dimension(:,:,:),
intent(inout),
optional :: data3d
473 integer :: variable_id
474 real(kind=default_precision),
dimension(:,:,:),
allocatable :: sdata
476 call check_status(nf90_inq_varid(ncid, key, variable_id))
478 if (.not.
present(data1d) .and. .not.
present(data3d))
return
480 if (
present(data1d))
then
481 call check_status(nf90_get_var(ncid, variable_id, data1d))
484 allocate(sdata(
size(data3d,1),
size(data3d,3),
size(data3d,2)))
485 call check_status(nf90_get_var(ncid, variable_id, sdata))
486 data3d(:,:,:)=reshape(sdata(:,:,:),(/
size(data3d,1),
size(data3d,2),
size(data3d,3)/))
495 integer,
intent(in) :: status
497 if (status /= nf90_noerr)
then
498 call log_log(log_error,
"NetCDF returned error code of "//trim(nf90_strerror(status)))