source: palm/trunk/SOURCE/parin.f90 @ 3241

Last change on this file since 3241 was 3240, checked in by Giersch, 6 years ago

Bufgix for an error message in case of restarts and user-defined profiles

  • Property svn:keywords set to Id
File size: 47.1 KB
Line 
1!> @file parin.f90
2!------------------------------------------------------------------------------!
3! This file is part of the PALM model system.
4!
5! PALM is free software: you can redistribute it and/or modify it under the
6! terms of the GNU General Public License as published by the Free Software
7! Foundation, either version 3 of the License, or (at your option) any later
8! version.
9!
10! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
11! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
12! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
13!
14! You should have received a copy of the GNU General Public License along with
15! PALM. If not, see <http://www.gnu.org/licenses/>.
16!
17! Copyright 1997-2018 Leibniz Universitaet Hannover
18!------------------------------------------------------------------------------!
19!
20! Current revisions:
21! -----------------
22!
23!
24! Former revisions:
25! -----------------
26! $Id: parin.f90 3240 2018-09-12 12:04:40Z raasch $
27! A check that controls the number of user-defined profiles on the restart file
28! with the one given for the current run has been added.
29!
30! 3204 2018-08-23 10:49:05Z raasch
31! additional check for nz
32!
33! 3184 2018-07-27 17:13:09Z suehring
34! Bugfix, comment setting of chemistry and passive scalar top boundary condition
35! in case of offline nesting
36!
37! 3183 2018-07-27 14:25:55Z suehring
38! Rename variables and boundary conditions in mesoscale-offline nesting mode
39!
40! 3182 2018-07-27 13:36:03Z suehring
41! Added multi agent system
42!
43! 3157 2018-07-19 21:08:49Z maronga
44! added use_free_convection_scaling
45!
46! 3083 2018-06-19 14:03:12Z gronemeier
47! Added rans_const_c and rans_const_sigma as input parameters (TG)
48!
49! 3065 2018-06-12 07:03:02Z Giersch
50! New initialization parameters added
51!
52! 3049 2018-05-29 13:52:36Z Giersch
53! Error messages revised
54!
55! 3045 2018-05-28 07:55:41Z Giersch
56! z_max_do2d removed, error messages revised
57!
58! 2995 2018-04-19 12:13:16Z Giersch
59! time_since_reference_point must be calculated/initialized before the first 
60! call of functions related to the radiation model which occur in
61! time_integration_spinup or time_integration
62!
63! 2980 2018-04-17 15:19:27Z suehring
64! Revise message call
65!
66! 2975 2018-04-16 15:22:20Z suehring
67! - Informative message when initializing_actions has been changed
68!   to set_constant_profile in child domain
69! - Change location in message call
70!
71! 2967 2018-04-13 11:22:08Z raasch
72! bugfix: missing parallel cpp-directives added
73!
74! 2941 2018-04-03 11:54:58Z kanani
75! Fix for spinup in case of restart run
76!
77! 2938 2018-03-27 15:52:42Z suehring
78! Change initialization in case child domain should be initialized with Inifor.
79!
80! 2936 2018-03-27 14:49:27Z suehring
81! inipar renamed to initialization_parameters.
82! d3par renamed to runtime_parameters.
83!
84! 2921 2018-03-22 15:05:23Z Giersch
85! Activation of spinup has been moved from lsm/usm_parin to parin itself
86!
87! 2906 2018-03-19 08:56:40Z Giersch
88! ENVIRONMENT variables read/write_svf has been added
89!
90! 2894 2018-03-15 09:17:58Z Giersch
91! read_var_list has been renamed to rrd_global, all module related _parin
92! routines are called before reading the global restart data to overwrite them
93! in case of restart runs
94!
95! 2881 2018-03-13 16:24:40Z suehring
96! Added flag for switching on/off calculation of soil moisture
97!
98! 2849 2018-03-05 10:49:33Z Giersch
99! Position of d3par namelist in parameter file is unimportant now
100!
101! 2826 2018-02-21 12:39:28Z Giersch
102! Bugfix in setting the default boundary conditions for nest domains
103!
104! 2817 2018-02-19 16:32:21Z knoop
105! Preliminary gust module interface implemented
106!
107! 2773 2018-01-30 14:12:54Z suehring
108! Nesting for chemical species implemented
109!
110! 2766 2018-01-22 17:17:47Z kanani
111! Removed preprocessor directive __chem
112!
113! 2718 2018-01-02 08:49:38Z maronga
114! Corrected "Former revisions" section
115!
116! 2696 2017-12-14 17:12:51Z kanani
117! Change in file header (GPL part)
118! Implementation of uv exposure model (FK)
119! Added rans_mode and turbulence_closure to inipar (TG)
120! Implementation of chemistry module
121! Sorting of USE list (FK)
122! Forcing implemented, and initialization with inifor (MS)
123!
124! 2600 2017-11-01 14:11:20Z raasch
125! some comments added and variables renamed concerning r2599
126!
127! 2599 2017-11-01 13:18:45Z hellstea
128! The i/o grouping is updated to work correctly also in nested runs.
129!
130! 2575 2017-10-24 09:57:58Z maronga
131! Renamed phi -> latitude, added longitude
132!
133! 2563 2017-10-19 15:36:10Z Giersch
134! Changed position where restart files are closed.
135!
136! 2550 2017-10-16 17:12:01Z boeske
137! Added complex_terrain
138!
139! 2544 2017-10-13 18:09:32Z maronga
140! Moved day_of_year_init and time_utc_init to inipar.
141!
142! 2397 2017-09-04 16:22:48Z suehring
143! Enable initialization of 3d model by user in the child domain.
144!
145! 2375 2017-08-29 14:10:28Z schwenkel
146! Added aerosol initialization for bulk microphysics
147!
148! 2372 2017-08-25 12:37:32Z sward
149! y_shift added to namelist
150!
151! 2365 2017-08-21 14:59:59Z kanani
152! Vertical grid nesting: add vnest_start_time to d3par (SadiqHuq)
153!
154! 2339 2017-08-07 13:55:26Z gronemeier
155! corrected timestamp in header
156!
157! 2338 2017-08-07 12:15:38Z gronemeier
158! Modularize 1D model
159!
160! 2310 2017-07-11 09:37:02Z gronemeier
161! Bugfix: re-arranged call for error messages for ENVPAR file
162!
163! 2304 2017-07-04 14:35:55Z suehring
164! Bugfix, enable restarts for child domain.
165!
166! 2298 2017-06-29 09:28:18Z raasch
167! -return_addres, return_username in ENVPAR, -cross_ts_uymax, cross_ts_uymin in
168! d3par
169!
170! 2296 2017-06-28 07:53:56Z maronga
171! Added parameters for model spinup
172!
173! 2292 2017-06-20 09:51:42Z schwenkel
174! Implementation of new microphysic scheme: cloud_scheme = 'morrison'
175! includes two more prognostic equations for cloud drop concentration (nc) 
176! and cloud water content (qc).
177!
178! 2267 2017-06-09 09:33:25Z gronemeier
179! Bugfix: removed skipping of reading namelists in case of omitted d3par
180!
181! 2259 2017-06-08 09:09:11Z gronemeier
182! Implemented synthetic turbulence generator
183!
184! 2233 2017-05-30 18:08:54Z suehring
185!
186! 2232 2017-05-30 17:47:52Z suehring
187! typo corrected
188! +wall_salinityflux
189! +tunnel_height, tunnel_lenght, tunnel_width_x, tunnel_width_y,
190!  tunnel_wall_depth
191!
192! 2118 2017-01-17 16:38:49Z raasch
193! -background_communication from inipar
194!
195! 2050 2016-11-08 15:00:55Z gronemeier
196! Implement turbulent outflow condition
197!
198! 2037 2016-10-26 11:15:40Z knoop
199! Anelastic approximation implemented
200!
201! 2035 2016-10-24 15:06:17Z suehring
202! Remove check for npex and npey in nesting case
203!
204! 2011 2016-09-19 17:29:57Z kanani
205! Added flag lsf_exception to allow explicit enabling of large scale forcing
206! together with buildings on flat terrain.
207!
208! 2007 2016-08-24 15:47:17Z kanani
209! Added call to urban surface model for reading of &urban_surface_par
210!
211! 2004 2016-08-24 10:25:59Z suehring
212! Humidity and passive scalar treated separately in nesting mode
213!
214! 2000 2016-08-20 18:09:15Z knoop
215! Forced header and separation lines into 80 columns
216!
217! 1992 2016-08-12 15:14:59Z suehring
218! +top_scalarflux
219!
220! 1960 2016-07-12 16:34:24Z suehring
221! Allocation of s_init
222!
223! 1957 2016-07-07 10:43:48Z suehring
224! flight module added
225!
226! 1955 2016-07-01 12:38:59Z hellstea
227! The parameter intializating_actions is set to 'set_constant_profiles for
228! all nest domains in order to make sure that diagnostic variables are properly
229! initialized for nest domains. Prognostic variables are later initialized by
230! interpolation from the parent domain.
231!
232! 1917 2016-05-27 14:28:12Z witha
233! Initial version of purely vertical nesting introduced.
234!
235! 1914 2016-05-26 14:44:07Z witha
236! Added call to wind turbine model for reading of &wind_turbine_par
237!
238! 1849 2016-04-08 11:33:18Z hoffmann
239! Adapted for modularization of microphysics
240!
241! 1833 2016-04-07 14:23:03Z raasch
242! call of spectra_parin
243!
244! 1831 2016-04-07 13:15:51Z hoffmann
245! turbulence renamed collision_turbulence, drizzle renamed
246! cloud_water_sedimentation
247! curvature_solution_effects removed
248!
249! 1826 2016-04-07 12:01:39Z maronga
250! Added call to radiation model for reading of &radiation_par.
251! Added call to plant canopy model for reading of &canopy_par.
252!
253! 1817 2016-04-06 15:44:20Z maronga
254! Added call to land surface model for reading of &lsm_par
255!
256! 1804 2016-04-05 16:30:18Z maronga
257! Removed code for parameter file check (__check)
258!
259! 1783 2016-03-06 18:36:17Z raasch
260! +netcdf_deflate in d3par, netcdf module and variable names changed
261!
262! 1764 2016-02-28 12:45:19Z raasch
263! cpp-statements for nesting removed, explicit settings of boundary conditions
264! in nest domains,
265! bugfix: npex/npey message moved from inipar to d3par
266! bugfix: check of lateral boundary conditions from check_parameters to here,
267! because they will be already used in init_pegrid and init_grid
268!
269! 1762 2016-02-25 12:31:13Z hellstea
270! Introduction of nested domain feature
271!
272! 1691 2015-10-26 16:17:44Z maronga
273! Added parameter most_method. Renamed prandtl_layer to constant_flux_layer.
274!
275! 1682 2015-10-07 23:56:08Z knoop
276! Code annotations made doxygen readable
277!
278! 1560 2015-03-06 10:48:54Z keck
279! +recycling_yshift
280!
281! 1496 2014-12-02 17:25:50Z maronga
282! Renamed: "radiation -> "cloud_top_radiation"
283!
284! 1484 2014-10-21 10:53:05Z kanani
285! Changes due to new module structure of the plant canopy model:
286!   canopy-model related parameters moved to new package canopy_par in
287!   subroutine package_parin
288!
289! 1429 2014-07-15 12:53:45Z knoop
290! +ensemble_member_nr to prepare the random_generator for ensemble runs
291!
292! 1402 2014-05-09 14:25:13Z raasch
293! location messages modified, batch_job included in envpar-NAMELIST
294!
295! 1384 2014-05-02 14:31:06Z raasch
296! location messages added
297!
298! 1365 2014-04-22 15:03:56Z boeske
299! Usage of large scale forcing enabled:
300! +use_subsidence_tendencies
301!
302! 1361 2014-04-16 15:17:48Z hoffmann
303! +call_microphysics_at_all_substeps
304!
305! 1359 2014-04-11 17:15:14Z hoffmann
306! REAL constants provided with KIND-attribute
307!
308! 1353 2014-04-08 15:21:23Z heinze
309! REAL constants provided with KIND-attribute
310!
311! 1327 2014-03-21 11:00:16Z raasch
312! -data_output_format, do3d_compress, do3d_comp_prec
313!
314! 1320 2014-03-20 08:40:49Z raasch
315! ONLY-attribute added to USE-statements,
316! kind-parameters added to all INTEGER and REAL declaration statements,
317! kinds are defined in new module kinds,
318! old module precision_kind is removed,
319! revision history before 2012 removed,
320! comment fields (!:) to be used for variable explanations added to
321! all variable declaration statements
322!
323! 1318 2014-03-17 13:35:16Z raasch
324! +cpu_log_barrierwait in d3par
325!
326! 1301 2014-03-06 13:29:46Z heinze
327! +large_scale_subsidence
328!
329! 1241 2013-10-30 11:36:58Z heinze
330! +nudging
331! +large_scale_forcing
332!
333! 1216 2013-08-26 09:31:42Z raasch
334! +transpose_compute_overlap in inipar
335!
336! 1195 2013-07-01 12:27:57Z heinze
337! Bugfix: allocate ref_state
338!
339! 1179 2013-06-14 05:57:58Z raasch
340! +reference_state in inipar
341!
342! 1159 2013-05-21 11:58:22Z fricke
343! +use_cmax
344!
345! 1128 2013-04-12 06:19:32Z raasch
346! +background_communication in inipar
347!
348! 1115 2013-03-26 18:16:16Z hoffmann
349! unused variables removed
350!
351! 1092 2013-02-02 11:24:22Z raasch
352! unused variables removed
353!
354! 1065 2012-11-22 17:42:36Z hoffmann
355! +nc, c_sedimentation, limiter_sedimentation, turbulence
356! -mu_constant, mu_constant_value
357!
358! 1053 2012-11-13 17:11:03Z hoffmann
359! necessary expansions according to the two new prognostic equations (nr, qr)
360! of the two-moment cloud physics scheme and steering parameters:
361! +*_init, *_surface, *_surface_initial_change, *_vertical_gradient,
362! +*_vertical_gradient_level, surface_waterflux_*,
363! +cloud_scheme, drizzle, mu_constant, mu_constant_value, ventilation_effect
364!
365! 1036 2012-10-22 13:43:42Z raasch
366! code put under GPL (PALM 3.9)
367!
368! 1015 2012-09-27 09:23:24Z raasch
369! -adjust_mixing_length
370!
371! 1003 2012-09-14 14:35:53Z raasch
372! -grid_matching
373!
374! 1001 2012-09-13 14:08:46Z raasch
375! -cut_spline_overshoot, long_filter_factor, overshoot_limit_*, ups_limit_*
376!
377! 996 2012-09-07 10:41:47Z raasch
378! -use_prior_plot1d_parameters
379!
380! 978 2012-08-09 08:28:32Z fricke
381! -km_damp_max, outflow_damping_width
382! +pt_damping_factor, pt_damping_width
383! +z0h_factor
384!
385! 964 2012-07-26 09:14:24Z raasch
386! -cross_normalized_x, cross_normalized_y, cross_xtext, z_max_do1d,
387! z_max_do1d_normalized
388!
389! 940 2012-07-09 14:31:00Z raasch
390! +neutral in inipar
391!
392! 927 2012-06-06 19:15:04Z raasch
393! +masking_method in inipar
394!
395! 824 2012-02-17 09:09:57Z raasch
396! +curvature_solution_effects in inipar
397!
398! 809 2012-01-30 13:32:58Z maronga
399! Bugfix: replaced .AND. and .NOT. with && and ! in the preprocessor directives
400!
401! 807 2012-01-25 11:53:51Z maronga
402! New cpp directive "__check" implemented which is used by check_namelist_files
403!
404! Revision 1.1  1997/07/24 11:22:50  raasch
405! Initial revision
406!
407!
408! Description:
409! ------------
410!> This subroutine reads variables controling the run from the NAMELIST files
411!------------------------------------------------------------------------------!
412 SUBROUTINE parin
413 
414
415    USE arrays_3d,                                                             &
416        ONLY:  pt_init, q_init, ref_state, s_init, sa_init,                    &     
417               ug, u_init, v_init, vg
418
419    USE chemistry_model_mod,                                                   &
420        ONLY:  chem_parin
421       
422    USE chem_modules
423
424    USE control_parameters
425
426    USE cpulog,                                                                &
427        ONLY:  cpu_log_barrierwait
428
429    USE date_and_time_mod,                                                     &
430        ONLY:  day_of_year_init, time_utc_init
431
432    USE dvrp_variables,                                                        &
433        ONLY:  local_dvrserver_running
434
435    USE flight_mod,                                                            &
436        ONLY:  flight_parin
437
438    USE grid_variables,                                                        &
439        ONLY:  dx, dy
440
441    USE gust_mod,                                                              &
442        ONLY: gust_parin
443
444    USE indices,                                                               &
445        ONLY:  nx, ny, nz
446
447    USE kinds
448
449    USE land_surface_model_mod,                                                &
450        ONLY: lsm_parin
451
452    USE microphysics_mod,                                                      &
453        ONLY:  c_sedimentation, cloud_water_sedimentation,                     &
454               collision_turbulence, curvature_solution_effects_bulk,          &
455               dry_aerosol_radius, limiter_sedimentation, na_init, nc_const,   &
456               sigma_bulk, ventilation_effect
457
458    USE model_1d_mod,                                                          &
459        ONLY:  damp_level_1d, dt_pr_1d, dt_run_control_1d, end_time_1d
460
461    USE multi_agent_system_mod,                                                &
462        ONLY:  mas_parin
463
464    USE netcdf_interface,                                                      &
465        ONLY:  netcdf_data_format, netcdf_deflate, netcdf_precision
466
467    USE pegrid
468               
469    USE plant_canopy_model_mod,                                                &
470         ONLY: pcm_parin
471
472    USE pmc_interface,                                                         &
473        ONLY:  nested_run, nesting_mode
474
475    USE profil_parameter,                                                      &
476        ONLY:  cross_profiles, profile_columns, profile_rows
477
478    USE progress_bar,                                                          &
479        ONLY :  batch_job
480
481    USE radiation_model_mod,                                                   &
482        ONLY: radiation_parin
483
484    USE read_restart_data_mod,                                                 &
485        ONLY:  rrd_global     
486
487    USE spectra_mod,                                                           &
488        ONLY :  spectra_parin
489
490    USE statistics,                                                            &
491        ONLY:  hom, hom_sum, pr_palm, region, statistic_regions
492
493    USE synthetic_turbulence_generator_mod,                                    &
494        ONLY:  stg_parin
495
496    USE turbulence_closure_mod,                                                &
497        ONLY:  rans_const_c, rans_const_sigma
498
499    USE urban_surface_mod,                                                     &
500        ONLY: usm_parin
501
502    USE uv_exposure_model_mod,                                                 &
503        ONLY:  uvem_parin
504
505    USE vertical_nesting_mod,                                                  &
506        ONLY:  vnest_start_time
507
508    USE wind_turbine_model_mod,                                                &
509        ONLY:  wtm_parin
510
511
512    IMPLICIT NONE
513
514    CHARACTER (LEN=80) ::  line  !< dummy string that contains the current line of the parameter file
515
516    INTEGER(iwp) ::  global_id      !< process id with respect to MPI_COMM_WORLD
517    INTEGER(iwp) ::  global_procs   !< # of procs with respect to MPI_COMM_WORLD
518    INTEGER(iwp) ::  i              !<
519    INTEGER(iwp) ::  ioerr          !< error flag for open/read/write
520
521    NAMELIST /inipar/  aerosol_bulk, alpha_surface, approximation, bc_e_b,     &
522                       bc_lr, bc_ns, bc_p_b, bc_p_t, bc_pt_b, bc_pt_t, bc_q_b, &
523             bc_q_t,bc_s_b, bc_s_t, bc_sa_t, bc_uv_b, bc_uv_t,                 &
524             bottom_salinityflux, building_height, building_length_x,          &
525             building_length_y, building_wall_left, building_wall_south,       &
526             calc_soil_moisture_during_spinup,                                 &
527             call_psolver_at_all_substeps, call_microphysics_at_all_substeps,  &
528             canyon_height,                                                    &
529             canyon_width_x, canyon_width_y, canyon_wall_left,                 &
530             canyon_wall_south, c_sedimentation, cfl_factor, cloud_droplets,   &
531             cloud_physics, cloud_scheme, cloud_top_radiation,                 &
532             cloud_water_sedimentation,                                        &
533             collective_wait, collision_turbulence, complex_terrain,           &
534             conserve_volume_flow,                                             &
535             conserve_volume_flow_mode, constant_flux_layer,                   &
536             coupling_start_time, curvature_solution_effects_bulk,             &
537             cycle_mg, damp_level_1d,                                          &
538             data_output_during_spinup,                                        &
539             day_of_year_init,                                                 &
540             dissipation_1d,                                                   &
541             dp_external, dp_level_b, dp_smooth, dpdxy, dry_aerosol_radius,    &
542             dt, dt_pr_1d, dt_run_control_1d, dt_spinup, dx, dy, dz, dz_max,   &
543             dz_stretch_factor, dz_stretch_level, dz_stretch_level_start,      &
544             dz_stretch_level_end, end_time_1d, ensemble_member_nr, e_init,    &
545             e_min, fft_method, flux_input_mode, flux_output_mode,             &
546             galilei_transformation, humidity,                                 &
547             inflow_damping_height, inflow_damping_width,                      &
548             inflow_disturbance_begin, inflow_disturbance_end,                 &
549             initializing_actions, km_constant,                                &
550             large_scale_forcing, large_scale_subsidence, latitude,            &
551             limiter_sedimentation, longitude,                                 &
552             loop_optimization, lsf_exception, masking_method, mg_cycles,      &
553             mg_switch_to_pe0_level, mixing_length_1d, momentum_advec,         &
554             most_method, na_init, nc_const, nesting_offline,                  &
555             netcdf_precision, neutral, ngsrb,                                 &
556             nsor, nsor_ini, nudging, nx, ny, nz, ocean, omega, omega_sor,     &
557             outflow_source_plane, passive_scalar,                             &
558             prandtl_number, precipitation, psolver, pt_damping_factor,        &
559             pt_damping_width, pt_reference, pt_surface,                       &
560             pt_surface_initial_change, pt_vertical_gradient,                  &
561             pt_vertical_gradient_level, q_surface, q_surface_initial_change,  &
562             q_vertical_gradient, q_vertical_gradient_level,                   &
563             random_generator, random_heatflux, rans_const_c, rans_const_sigma,&
564             rans_mode,                                                        &
565             rayleigh_damping_factor, rayleigh_damping_height,                 &
566             recycling_width, recycling_yshift,                                &
567             reference_state, residual_limit,                                  &
568             roughness_length, sa_surface,                                     &
569             sa_vertical_gradient, sa_vertical_gradient_level, scalar_advec,   &
570             scalar_rayleigh_damping, sigma_bulk,                              &
571             spinup_time, spinup_pt_amplitude, spinup_pt_mean,                 &
572             statistic_regions, subs_vertical_gradient,                        &
573             subs_vertical_gradient_level, surface_heatflux, surface_pressure, &
574             surface_scalarflux, surface_waterflux,                            &
575             s_surface, s_surface_initial_change, s_vertical_gradient,         &
576             s_vertical_gradient_level, time_utc_init, timestep_scheme,        &
577             topography, topography_grid_convention, top_heatflux,             &
578             top_momentumflux_u, top_momentumflux_v, top_salinityflux,         &
579             top_scalarflux, transpose_compute_overlap,                        &
580             tunnel_height, tunnel_length, tunnel_width_x, tunnel_width_y,     &
581             tunnel_wall_depth, turbulence_closure,                            &
582             turbulent_inflow, turbulent_outflow,                              &
583             use_subsidence_tendencies, ug_surface, ug_vertical_gradient,      &
584             use_free_convection_scaling,                                      &
585             ug_vertical_gradient_level, use_surface_fluxes, use_cmax,         &
586             use_top_fluxes, use_ug_for_galilei_tr, use_upstream_for_tke,      &
587             uv_heights, u_bulk, u_profile, vg_surface, vg_vertical_gradient,  &
588             vg_vertical_gradient_level, v_bulk, v_profile, ventilation_effect,&
589             wall_adjustment, wall_heatflux, wall_humidityflux,                &
590             wall_salinityflux, wall_scalarflux, y_shift, zeta_max, zeta_min,  &
591             z0h_factor
592
593    NAMELIST /initialization_parameters/  aerosol_bulk, alpha_surface,         &
594             approximation, bc_e_b,                                            &
595             bc_lr, bc_ns, bc_p_b, bc_p_t, bc_pt_b, bc_pt_t, bc_q_b,           &
596             bc_q_t,bc_s_b, bc_s_t, bc_sa_t, bc_uv_b, bc_uv_t,                 &
597             bottom_salinityflux, building_height, building_length_x,          &
598             building_length_y, building_wall_left, building_wall_south,       &
599             calc_soil_moisture_during_spinup,                                 &
600             call_psolver_at_all_substeps, call_microphysics_at_all_substeps,  &
601             canyon_height,                                                    &
602             canyon_width_x, canyon_width_y, canyon_wall_left,                 &
603             canyon_wall_south, c_sedimentation, cfl_factor, cloud_droplets,   &
604             cloud_physics, cloud_scheme, cloud_top_radiation,                 &
605             cloud_water_sedimentation,                                        &
606             collective_wait, collision_turbulence, complex_terrain,           &
607             conserve_volume_flow,                                             &
608             conserve_volume_flow_mode, constant_flux_layer,                   &
609             coupling_start_time, curvature_solution_effects_bulk,             &
610             cycle_mg, damp_level_1d,                                          &
611             data_output_during_spinup,                                        &
612             day_of_year_init,                                                 &
613             dissipation_1d,                                                   &
614             dp_external, dp_level_b, dp_smooth, dpdxy, dry_aerosol_radius,    &
615             dt, dt_pr_1d, dt_run_control_1d, dt_spinup, dx, dy, dz, dz_max,   &
616             dz_stretch_factor, dz_stretch_level, dz_stretch_level_start,      &
617             dz_stretch_level_end, end_time_1d, ensemble_member_nr, e_init,    &
618             e_min, fft_method, flux_input_mode, flux_output_mode,             &
619             galilei_transformation, humidity,                                 &
620             inflow_damping_height, inflow_damping_width,                      &
621             inflow_disturbance_begin, inflow_disturbance_end,                 &
622             initializing_actions, km_constant,                                &
623             large_scale_forcing, large_scale_subsidence, latitude,            &
624             limiter_sedimentation, longitude,                                 &
625             loop_optimization, lsf_exception, masking_method, mg_cycles,      &
626             mg_switch_to_pe0_level, mixing_length_1d, momentum_advec,         &
627             most_method, na_init, nc_const, nesting_offline,                  &
628             netcdf_precision, neutral, ngsrb,                                 &
629             nsor, nsor_ini, nudging, nx, ny, nz, ocean, omega, omega_sor,     &
630             outflow_source_plane, passive_scalar,                             &
631             prandtl_number, precipitation, psolver, pt_damping_factor,        &
632             pt_damping_width, pt_reference, pt_surface,                       &
633             pt_surface_initial_change, pt_vertical_gradient,                  &
634             pt_vertical_gradient_level, q_surface, q_surface_initial_change,  &
635             q_vertical_gradient, q_vertical_gradient_level,                   &
636             random_generator, random_heatflux, rans_const_c, rans_const_sigma,&
637             rans_mode,                                                        &
638             rayleigh_damping_factor, rayleigh_damping_height,                 &
639             recycling_width, recycling_yshift,                                &
640             reference_state, residual_limit,                                  &
641             roughness_length, sa_surface,                                     &
642             sa_vertical_gradient, sa_vertical_gradient_level, scalar_advec,   &
643             scalar_rayleigh_damping, sigma_bulk,                              &
644             spinup_time, spinup_pt_amplitude, spinup_pt_mean,                 &
645             statistic_regions, subs_vertical_gradient,                        &
646             subs_vertical_gradient_level, surface_heatflux, surface_pressure, &
647             surface_scalarflux, surface_waterflux,                            &
648             s_surface, s_surface_initial_change, s_vertical_gradient,         &
649             s_vertical_gradient_level, time_utc_init, timestep_scheme,        &
650             topography, topography_grid_convention, top_heatflux,             &
651             top_momentumflux_u, top_momentumflux_v, top_salinityflux,         &
652             top_scalarflux, transpose_compute_overlap,                        &
653             tunnel_height, tunnel_length, tunnel_width_x, tunnel_width_y,     &
654             tunnel_wall_depth, turbulence_closure,                            &
655             turbulent_inflow, turbulent_outflow,                              &
656             use_subsidence_tendencies, ug_surface, ug_vertical_gradient,      &
657             ug_vertical_gradient_level, use_surface_fluxes, use_cmax,         &
658             use_top_fluxes, use_ug_for_galilei_tr, use_upstream_for_tke,      &
659             use_free_convection_scaling,                                      &
660             uv_heights, u_bulk, u_profile, vg_surface, vg_vertical_gradient,  &
661             vg_vertical_gradient_level, v_bulk, v_profile, ventilation_effect,&
662             wall_adjustment, wall_heatflux, wall_humidityflux,                &
663             wall_salinityflux, wall_scalarflux, y_shift, zeta_max, zeta_min,  &
664             z0h_factor
665             
666    NAMELIST /d3par/  averaging_interval, averaging_interval_pr,               &
667             cpu_log_barrierwait, create_disturbances,                         &
668             cross_profiles, data_output, data_output_masks,                   &
669             data_output_pr, data_output_2d_on_each_pe, disturbance_amplitude, &
670             disturbance_energy_limit, disturbance_level_b,                    &
671             disturbance_level_t, do2d_at_begin, do3d_at_begin,                &
672             dt, dt_averaging_input, dt_averaging_input_pr,                    &
673             dt_coupling, dt_data_output, dt_data_output_av, dt_disturb,       &
674             dt_domask, dt_dopr, dt_dopr_listing, dt_dots, dt_do2d_xy,         &
675             dt_do2d_xz, dt_do2d_yz, dt_do3d, dt_max, dt_restart,              &
676             dt_run_control,end_time, force_print_header, mask_scale_x,        &
677             mask_scale_y, mask_scale_z, mask_x, mask_y, mask_z, mask_x_loop,  &
678             mask_y_loop, mask_z_loop, netcdf_data_format, netcdf_deflate,     &
679             normalizing_region, npex, npey, nz_do3d,                          &
680             precipitation_amount_interval, profile_columns, profile_rows,     &
681             restart_time, section_xy, section_xz, section_yz,                 &
682             skip_time_data_output, skip_time_data_output_av, skip_time_dopr,  &
683             skip_time_do2d_xy, skip_time_do2d_xz, skip_time_do2d_yz,          &
684             skip_time_do3d, skip_time_domask, synchronous_exchange,           &
685             termination_time_needed, vnest_start_time
686
687    NAMELIST /runtime_parameters/  averaging_interval, averaging_interval_pr,  &
688             cpu_log_barrierwait, create_disturbances,                         &
689             cross_profiles, data_output, data_output_masks,                   &
690             data_output_pr, data_output_2d_on_each_pe, disturbance_amplitude, &
691             disturbance_energy_limit, disturbance_level_b,                    &
692             disturbance_level_t, do2d_at_begin, do3d_at_begin,                &
693             dt, dt_averaging_input, dt_averaging_input_pr,                    &
694             dt_coupling, dt_data_output, dt_data_output_av, dt_disturb,       &
695             dt_domask, dt_dopr, dt_dopr_listing, dt_dots, dt_do2d_xy,         &
696             dt_do2d_xz, dt_do2d_yz, dt_do3d, dt_max, dt_restart,              &
697             dt_run_control,end_time, force_print_header, mask_scale_x,        &
698             mask_scale_y, mask_scale_z, mask_x, mask_y, mask_z, mask_x_loop,  &
699             mask_y_loop, mask_z_loop, netcdf_data_format, netcdf_deflate,     &
700             normalizing_region, npex, npey, nz_do3d,                          &
701             precipitation_amount_interval, profile_columns, profile_rows,     &
702             restart_time, section_xy, section_xz, section_yz,                 &
703             skip_time_data_output, skip_time_data_output_av, skip_time_dopr,  &
704             skip_time_do2d_xy, skip_time_do2d_xz, skip_time_do2d_yz,          &
705             skip_time_do3d, skip_time_domask, synchronous_exchange,           &
706             termination_time_needed, vnest_start_time
707
708    NAMELIST /envpar/  batch_job, host, local_dvrserver_running,               &
709                       maximum_cpu_time_allowed, maximum_parallel_io_streams,  &
710                       read_svf, revision, run_identifier, tasks_per_node,     &
711                       write_binary, write_svf
712
713!
714!-- First read values of environment variables (this NAMELIST file is
715!-- generated by palmrun)
716    CALL location_message( 'reading environment parameters from ENVPAR', .FALSE. )
717
718    OPEN ( 90, FILE='ENVPAR', STATUS='OLD', FORM='FORMATTED', IOSTAT=ioerr )
719
720    IF ( ioerr /= 0 )  THEN
721       message_string = 'local file ENVPAR not found' //                       &
722                        '&some variables for steering may not be properly set'
723       CALL message( 'parin', 'PA0276', 0, 1, 0, 6, 0 )
724    ELSE
725       READ ( 90, envpar, IOSTAT=ioerr )
726       IF ( ioerr < 0 )  THEN
727          message_string = 'no envpar-NAMELIST found in local file '  //       &
728                           'ENVPAR& or some variables for steering may '  //   &
729                           'not be properly set'
730          CALL message( 'parin', 'PA0278', 0, 1, 0, 6, 0 )
731       ELSEIF ( ioerr > 0 )  THEN
732          message_string = 'errors in local file ENVPAR' //                    &
733                           '&some variables for steering may not be properly set'
734          CALL message( 'parin', 'PA0277', 0, 1, 0, 6, 0 )
735       ENDIF
736       CLOSE ( 90 )
737    ENDIF
738
739    CALL location_message( 'finished', .TRUE. )
740!
741!-- Calculate the number of groups into which parallel I/O is split.
742!-- The default for files which are opened by all PEs (or where each
743!-- PE opens his own independent file) is, that all PEs are doing input/output
744!-- in parallel at the same time. This might cause performance or even more
745!-- severe problems depending on the configuration of the underlying file
746!-- system.
747!-- Calculation of the number of blocks and the I/O group must be based on all
748!-- PEs involved in this run. Since myid and numprocs are related to the
749!-- comm2d communicator, which gives only a subset of all PEs in case of
750!-- nested runs, that information must be inquired again from the global
751!-- communicator.
752!-- First, set the default:
753#if defined( __parallel )
754    CALL MPI_COMM_RANK( MPI_COMM_WORLD, global_id, ierr )
755    CALL MPI_COMM_SIZE( MPI_COMM_WORLD, global_procs, ierr )
756#else
757    global_id    = 0
758    global_procs = 1
759#endif
760    IF ( maximum_parallel_io_streams == -1  .OR.                               &
761         maximum_parallel_io_streams > global_procs )  THEN
762       maximum_parallel_io_streams = global_procs
763    ENDIF
764!
765!-- Now calculate the number of io_blocks and the io_group to which the
766!-- respective PE belongs. I/O of the groups is done in serial, but in parallel
767!-- for all PEs belonging to the same group.
768    io_blocks = global_procs / maximum_parallel_io_streams
769    io_group  = MOD( global_id+1, io_blocks )
770   
771    CALL location_message( 'reading NAMELIST parameters from PARIN', .FALSE. )
772!
773!-- Data is read in parallel by groups of PEs
774    DO  i = 0, io_blocks-1
775       IF ( i == io_group )  THEN
776
777!
778!--       Open the NAMELIST-file which is send with this job
779          CALL check_open( 11 )
780
781!
782!--       Read the control parameters for initialization.
783!--       The namelist "inipar" must be provided in the NAMELIST-file.
784          READ ( 11, initialization_parameters, ERR=10, END=11 )
785
786          GOTO 12
787
788 10       message_string = 'errors in initialization_parameters & or no ' //  &
789                           'initialization_parameters-namelist ' //           &
790                           'found (CRAY-machines only)'
791          CALL message( 'parin', 'PA0271', 1, 2, 0, 6, 0 )
792
793 11       REWIND ( 11 )
794          READ ( 11, inipar, ERR=13, END=14 )
795 
796          message_string = 'namelist inipar is deprecated and will be ' //    &
797                          'removed in near future. & Please use namelist ' // &
798                          'initialization_parameters instead'
799          CALL message( 'parin', 'PA0017', 0, 1, 0, 6, 0 )
800 
801          GOTO 12
802 
803 13       message_string = 'errors in inipar & or no inipar-namelist ' //      &
804                           'found (CRAY-machines only)'
805          CALL message( 'parin', 'PA0271', 1, 2, 0, 6, 0 )
806         
807 14       message_string = 'no initialization_parameters-namelist found'
808          CALL message( 'parin', 'PA0272', 1, 2, 0, 6, 0 )
809
810!
811!--       Try to read runtime parameters given by the user for this run
812!--       (namelist "runtime_parameters"). The namelist "runtime_parmeters"   
813!--       can be omitted. In that case default values are used for the         
814!--       parameters.
815 12       line = ' '
816
817          REWIND ( 11 )
818          line = ' '
819          DO   WHILE ( INDEX( line, '&runtime_parameters' ) == 0 )
820             READ ( 11, '(A)', END=20 )  line
821          ENDDO
822          BACKSPACE ( 11 )
823
824!
825!--       Read namelist
826          READ ( 11, runtime_parameters )
827
828          GOTO 21
829         
830 20       REWIND ( 11 )
831          line = ' '
832          DO   WHILE ( INDEX( line, '&d3par' ) == 0 )
833             READ ( 11, '(A)', END=21 )  line
834          ENDDO
835          BACKSPACE ( 11 )
836 
837 !
838!--       Read namelist
839          READ ( 11, d3par )
840 
841          message_string = 'namelist d3par is deprecated and will be ' //      &
842                          'removed in near future. &Please use namelist ' //   &
843                          'runtime_parameters instead'
844          CALL message( 'parin', 'PA0487', 0, 1, 0, 6, 0 )
845         
846 21       CONTINUE
847
848!
849!--       Check if land surface model is used and read &lsm_par if required
850          CALL lsm_parin
851
852!
853!--       Check if urban surface model is used and read &urban_surface_par if required
854          CALL usm_parin
855
856!
857!--       Check if spectra shall be calculated and read spectra_par if required
858          CALL spectra_parin
859
860!
861!--       Check if radiation model is used and read &radiation_par if required
862          CALL radiation_parin
863
864!
865!--       Check if gust module is used and read &gust_par if required
866          CALL gust_parin
867
868!
869!--       Check if multi agent system is used and read &agents_par if required
870          CALL mas_parin
871
872!
873!--       Check if plant canopy model is used and read &canopy_par if required
874          CALL pcm_parin
875
876!
877!--       Read control parameters for optionally used model software packages
878          CALL package_parin
879
880!
881!--       Check if wind turbine model is used and read &wind_turbine_par if
882!--       required
883          CALL wtm_parin
884!
885!--       Check if virtual flights should be carried out and read &flight_par
886!--       if required
887          CALL flight_parin
888!
889!--       Check if synthetic turbulence generator is used and read stg_par if
890!--       required
891          CALL stg_parin
892!
893!--       Read chemistry variables
894          CALL chem_parin
895!
896!--       Check if uv exposure model is used and read &uvexposure_par
897          CALL uvem_parin
898!
899!--       Read user-defined variables
900          CALL user_parin
901
902!
903!--       If required, read control parameters from restart file (produced by
904!--       a prior run). All PEs are reading from file created by PE0 (see
905!--       check_open)
906          IF ( TRIM( initializing_actions ) == 'read_restart_data' )  THEN
907
908             CALL rrd_global
909!
910!--          Increment the run count
911             runnr = runnr + 1
912!
913!--          In case of a restart run, the number of user-defined profiles on
914!--          the restart file (already stored in max_pr_user) has to match the
915!--          one given for the current run. max_pr_user_tmp is calculated in
916!--          user_parin and max_pr_user is read in via rrd_global.
917             IF ( max_pr_user /= max_pr_user_tmp )  THEN
918                WRITE( message_string, * ) 'the number of user-defined ',      &
919                      'profiles given in data_output_pr (', max_pr_user_tmp,   &
920                      ') does not match the one ',                             &
921                      'found in the restart file (', max_pr_user, ')'
922                CALL message( 'user_parin', 'UI0009', 1, 2, 0, 6, 0 )
923             ENDIF
924          ELSE
925             max_pr_user = max_pr_user_tmp
926          ENDIF
927
928!
929!--       Activate spinup
930          IF ( land_surface .OR. urban_surface )  THEN
931             IF ( spinup_time > 0.0_wp )  THEN
932                coupling_start_time = spinup_time
933                time_since_reference_point = simulated_time - coupling_start_time
934                IF ( spinup_pt_mean == 9999999.9_wp )  THEN
935                   spinup_pt_mean = pt_surface
936                ENDIF
937                end_time = end_time + spinup_time
938                IF ( TRIM( initializing_actions ) /= 'read_restart_data' )     &
939                   spinup = .TRUE.
940             ENDIF
941          ENDIF
942
943!
944!--       In case of nested runs, explicitly set nesting boundary conditions.
945!--       This will overwrite the user settings and basic defaults.
946!--       bc_lr and bc_ns always need to be cyclic for vertical nesting.
947          IF ( nested_run )  THEN
948             IF ( nesting_mode == 'vertical' )  THEN
949                IF (bc_lr /= 'cyclic' .OR. bc_ns /= 'cyclic' )  THEN
950                   WRITE ( message_string, *) 'bc_lr and bc_ns were set to ,', &
951                        'cyclic for vertical nesting'
952                   CALL message( 'parin', 'PA0428', 0, 0, 0, 6, 0 )
953                   bc_lr   = 'cyclic'
954                   bc_ns   = 'cyclic'
955                ENDIF
956                IF ( child_domain )  THEN
957                   bc_uv_t  = 'nested'
958                   bc_pt_t  = 'nested'
959                   bc_q_t   = 'nested'
960                   bc_s_t   = 'nested'
961                   bc_cs_t  = 'nested'
962                   bc_p_t   = 'neumann' 
963                ENDIF
964!
965!--          For other nesting modes only set boundary conditions for
966!--          nested domains.
967             ELSE
968                IF ( child_domain )  THEN
969                   bc_lr    = 'nested'
970                   bc_ns    = 'nested'
971                   bc_uv_t  = 'nested'
972                   bc_pt_t  = 'nested'
973                   bc_q_t   = 'nested'
974                   bc_s_t   = 'nested'
975                   bc_cs_t  = 'nested'
976                   bc_p_t   = 'neumann'
977                ENDIF
978             ENDIF
979          ENDIF
980!
981!--       Set boundary conditions also in case the model is offline-nested in
982!--       larger-scale models.
983          IF ( nesting_offline )  THEN
984             bc_lr    = 'nesting_offline'
985             bc_ns    = 'nesting_offline'
986             bc_uv_t  = 'nesting_offline'
987             bc_pt_t  = 'nesting_offline'
988             bc_q_t   = 'nesting_offline'
989           !  bc_s_t   = 'nesting_offline'  ! scalar boundary condition is not clear
990           !  bc_cs_t  = 'nesting_offline'  ! same for chemical species
991             bc_p_t   = 'neumann'
992          ENDIF
993
994!         
995!--       In case of nested runs, make sure that initializing_actions =
996!--       'set_constant_profiles' even though the constant-profiles
997!--       initializations for the prognostic variables will be overwritten
998!--       by pmci_child_initialize and pmci_parent_initialize. This is,
999!--       however, important e.g. to make sure that diagnostic variables
1000!--       are set properly. An exception is made in case of restart runs and
1001!--       if user decides to do everything by its own.
1002          IF ( child_domain  .AND.  .NOT. (                                    &
1003               TRIM( initializing_actions ) == 'read_restart_data'      .OR.   &
1004               TRIM( initializing_actions ) == 'set_constant_profiles'  .OR.   &
1005               TRIM( initializing_actions ) == 'by_user' ) )  THEN
1006             message_string = 'initializing_actions = ' //                     &
1007                              TRIM( initializing_actions ) // ' has been ' //  &
1008                              'changed to set_constant_profiles in child ' //  &
1009                              'domain.' 
1010             CALL message( 'parin', 'PA0492', 0, 0, 0, 6, 0 )
1011
1012             initializing_actions = 'set_constant_profiles'
1013          ENDIF           
1014!
1015!--       Check validity of lateral boundary conditions. This has to be done
1016!--       here because they are already used in init_pegrid and init_grid and
1017!--       therefore cannot be check in check_parameters
1018          IF ( bc_lr /= 'cyclic'  .AND.  bc_lr /= 'dirichlet/radiation'  .AND. &
1019               bc_lr /= 'radiation/dirichlet'  .AND.  bc_lr /= 'nested'  .AND. &
1020               bc_lr /= 'nesting_offline' )  THEN
1021             message_string = 'unknown boundary condition: bc_lr = "' // &
1022                              TRIM( bc_lr ) // '"'
1023             CALL message( 'parin', 'PA0049', 1, 2, 0, 6, 0 )
1024          ENDIF
1025          IF ( bc_ns /= 'cyclic'  .AND.  bc_ns /= 'dirichlet/radiation'  .AND. &
1026               bc_ns /= 'radiation/dirichlet'  .AND.  bc_ns /= 'nested'  .AND. &
1027               bc_ns /= 'nesting_offline' )  THEN
1028             message_string = 'unknown boundary condition: bc_ns = "' // &
1029                              TRIM( bc_ns ) // '"'
1030             CALL message( 'parin', 'PA0050', 1, 2, 0, 6, 0 )
1031          ENDIF
1032!
1033!--       Set internal variables used for speed optimization in if clauses
1034          IF ( bc_lr /= 'cyclic' )               bc_lr_cyc    = .FALSE.
1035          IF ( bc_lr == 'dirichlet/radiation' )  bc_lr_dirrad = .TRUE.
1036          IF ( bc_lr == 'radiation/dirichlet' )  bc_lr_raddir = .TRUE.
1037          IF ( bc_ns /= 'cyclic' )               bc_ns_cyc    = .FALSE.
1038          IF ( bc_ns == 'dirichlet/radiation' )  bc_ns_dirrad = .TRUE.
1039          IF ( bc_ns == 'radiation/dirichlet' )  bc_ns_raddir = .TRUE.
1040
1041!
1042!--       Definition of names of areas used for computing statistics. They must
1043!--       be defined at this place, because they are allowed to be redefined by
1044!--       the user in user_parin.
1045          region = 'total domain'
1046
1047!
1048!--       Check in case of initial run, if the grid point numbers are well
1049!--       defined and allocate some arrays which are already needed in
1050!--       init_pegrid or check_parameters. During restart jobs, these arrays
1051!--       will be allocated in rrd_global. All other arrays are allocated
1052!--       in init_3d_model.
1053          IF ( TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
1054
1055             IF ( nx <= 0 )  THEN
1056                WRITE( message_string, * ) 'no value or wrong value given',    &
1057                                           ' for nx: nx=', nx
1058                CALL message( 'parin', 'PA0273', 1, 2, 0, 6, 0 )
1059             ENDIF
1060             IF ( ny <= 0 )  THEN
1061                WRITE( message_string, * ) 'no value or wrong value given',    &
1062                                           ' for ny: ny=', ny
1063                CALL message( 'parin', 'PA0274', 1, 2, 0, 6, 0 )
1064             ENDIF
1065             IF ( nz <= 0 )  THEN
1066                WRITE( message_string, * ) 'no value or wrong value given',    &
1067                                           ' for nz: nz=', nz
1068                CALL message( 'parin', 'PA0275', 1, 2, 0, 6, 0 )
1069             ENDIF
1070
1071!
1072!--          As a side condition, routine flow_statistics require at least 14
1073!--          vertical grid levels (they are used to store time-series data)
1074!>           @todo   Remove this restriction
1075             IF ( nz < 14 )  THEN
1076                WRITE( message_string, * ) 'nz >= 14 is required'
1077                CALL message( 'parin', 'PA0362', 1, 2, 0, 6, 0 )
1078             ENDIF
1079
1080!
1081!--          ATTENTION: in case of changes to the following statement please
1082!--                  also check the allocate statement in routine rrd_global
1083             ALLOCATE( pt_init(0:nz+1), q_init(0:nz+1), s_init(0:nz+1),        &
1084                       ref_state(0:nz+1), sa_init(0:nz+1), ug(0:nz+1),         &
1085                       u_init(0:nz+1), v_init(0:nz+1), vg(0:nz+1),             &
1086                       hom(0:nz+1,2,pr_palm+max_pr_user,0:statistic_regions),  &
1087                       hom_sum(0:nz+1,pr_palm+max_pr_user,0:statistic_regions) )
1088
1089             hom = 0.0_wp
1090
1091          ENDIF
1092
1093!
1094!--       NAMELIST-file is not needed anymore
1095          CALL close_file( 11 )
1096
1097       ENDIF
1098#if defined( __parallel )
1099       CALL MPI_BARRIER( MPI_COMM_WORLD, ierr )
1100#endif
1101    ENDDO
1102
1103    CALL location_message( 'finished', .TRUE. )
1104
1105 END SUBROUTINE parin
Note: See TracBrowser for help on using the repository browser.