Skip to content

Commit 744ab4f

Browse files
authored
Make initial iterations compatible with boundary planes (Exawind#1651)
* enable initial iterations to be compatible with boundary planes * set up reg test for this capability * formatting
1 parent 412f015 commit 744ab4f

File tree

6 files changed

+88
-6
lines changed

6 files changed

+88
-6
lines changed

amr-wind/equation_systems/PDE.H

Lines changed: 3 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -137,7 +137,9 @@ public:
137137
{
138138
BL_PROFILE(
139139
"amr-wind::" + this->identifier() + "::pre_advection_actions");
140-
m_adv_op->preadvect(fstate, m_time.delta_t(), m_time.current_time());
140+
m_adv_op->preadvect(
141+
fstate, m_time.delta_t(),
142+
0.5 * (m_time.current_time() + m_time.new_time()));
141143
}
142144

143145
void compute_predictor_rhs(const DiffusionType difftype) override

amr-wind/equation_systems/PDEBase.cpp

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -94,7 +94,7 @@ void PDEMgr::prepare_boundaries()
9494
{
9595
// If state variables exist at NPH, fill their boundary cells
9696
const auto nph_time =
97-
m_sim.time().current_time() + 0.5 * m_sim.time().delta_t();
97+
0.5 * (m_sim.time().current_time() + m_sim.time().new_time());
9898
if (m_constant_density &&
9999
m_sim.repo().field_exists("density", FieldState::NPH)) {
100100
auto& nph_field =

amr-wind/equation_systems/icns/icns_advection.H

Lines changed: 2 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -208,7 +208,7 @@ struct AdvectionOp<ICNS, fvm::Godunov>
208208
}
209209

210210
// Populate boundaries (valid cells) using velocity BCs
211-
m_macproj_op.set_inflow_velocity(time + 0.5 * dt);
211+
m_macproj_op.set_inflow_velocity(time);
212212

213213
// MAC projection
214214
m_macproj_op(fstate, dt);
@@ -217,8 +217,7 @@ struct AdvectionOp<ICNS, fvm::Godunov>
217217
if (fvm::Godunov::nghost_state > 0) {
218218
amrex::Array<Field*, AMREX_SPACEDIM> mac_vel = {
219219
AMREX_D_DECL(&u_mac, &v_mac, &w_mac)};
220-
dof_field.fillpatch_sibling_fields(
221-
time + 0.5 * dt, u_mac.num_grow(), mac_vel);
220+
dof_field.fillpatch_sibling_fields(time, u_mac.num_grow(), mac_vel);
222221
}
223222

224223
for (int lev = 0; lev < repo.num_active_levels(); ++lev) {

test/CMakeLists.txt

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -403,6 +403,7 @@ endif()
403403

404404
if(AMR_WIND_ENABLE_NETCDF)
405405
add_test_red(abl_bndry_input abl_bndry_output)
406+
add_test_red(abl_bndry_input_init abl_bndry_output)
406407
add_test_red(abl_bndry_input_amr abl_bndry_output)
407408
add_test_red(abl_bndry_input_amr_inflow abl_bndry_output_amr_inflow)
408409
add_test_red(abl_bndry_input_amr_upper abl_bndry_output_amr_upper)
Lines changed: 80 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,80 @@
1+
#¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨#
2+
# SIMULATION STOP #
3+
#.......................................#
4+
time.stop_time = 22000.0 # Max (simulated) time to evolve
5+
time.max_step = 10 # Max number of time steps
6+
#¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨#
7+
# TIME STEP COMPUTATION #
8+
#.......................................#
9+
time.fixed_dt = 0.4 # Use this constant dt if > 0
10+
time.cfl = 0.95 # CFL factor
11+
#¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨#
12+
# INPUT AND OUTPUT #
13+
#.......................................#
14+
time.plot_interval = 10 # Steps between plot files
15+
time.checkpoint_interval = -1 # Steps between checkpoint files
16+
#¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨#
17+
# PHYSICS #
18+
#.......................................#
19+
incflo.gravity = 0. 0. -9.81 # Gravitational force (3D)
20+
incflo.density = 1.0 # Reference density
21+
incflo.use_godunov = 1
22+
incflo.diffusion_type = 2
23+
transport.viscosity = 1.0e-5
24+
transport.laminar_prandtl = 0.7
25+
transport.turbulent_prandtl = 0.3333
26+
transport.reference_temperature = 290.0
27+
turbulence.model = Smagorinsky
28+
Smagorinsky_coeffs.Cs = 0.135
29+
incflo.physics = ABL
30+
ICNS.source_terms = CoriolisForcing GeostrophicForcing
31+
CoriolisForcing.east_vector = 1.0 0.0 0.0
32+
CoriolisForcing.north_vector = 0.0 1.0 0.0
33+
CoriolisForcing.latitude = 90.0
34+
CoriolisForcing.rotational_time_period = 125663.706143592
35+
GeostrophicForcing.geostrophic_wind = 10.0 0.0 0.0
36+
incflo.velocity = 10.0 0.0 0.0
37+
ABL.temperature_heights = 0.0 2000.0
38+
ABL.temperature_values = 290.0 290.0
39+
ABL.perturb_temperature = false
40+
ABL.cutoff_height = 50.0
41+
ABL.perturb_velocity = true
42+
ABL.perturb_ref_height = 50.0
43+
ABL.Uperiods = 4.0
44+
ABL.Vperiods = 4.0
45+
ABL.deltaU = 1.0
46+
ABL.deltaV = 1.0
47+
ABL.kappa = .41
48+
ABL.surface_roughness_z0 = 0.01
49+
ABL.bndry_file = "../abl_bndry_output/bndry_file.nc"
50+
ABL.bndry_io_mode = 1
51+
ABL.bndry_var_names = velocity temperature
52+
ABL.bndry_output_format = netcdf
53+
#¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨#
54+
# ADAPTIVE MESH REFINEMENT #
55+
#.......................................#
56+
amr.n_cell = 48 48 48 # Grid cells at coarsest AMRlevel
57+
amr.max_level = 0 # Max AMR level in hierarchy
58+
#¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨#
59+
# GEOMETRY #
60+
#.......................................#
61+
geometry.prob_lo = 0. 0. 0. # Lo corner coordinates
62+
geometry.prob_hi = 1000. 1000. 1000. # Hi corner coordinates
63+
geometry.is_periodic = 0 0 0 # Periodicity x y z (0/1)
64+
# Boundary conditions
65+
xlo.type = "mass_inflow"
66+
xlo.density = 1.0
67+
xlo.temperature = 0.0
68+
xhi.type = "pressure_outflow"
69+
ylo.type = "mass_inflow"
70+
ylo.density = 1.0
71+
ylo.temperature = 0.0
72+
yhi.type = "pressure_outflow"
73+
zlo.type = "wall_model"
74+
zhi.type = "slip_wall"
75+
zhi.temperature_type = "fixed_gradient"
76+
zhi.temperature = 0.0
77+
#¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨#
78+
# VERBOSITY #
79+
#.......................................#
80+
incflo.verbose = 0 # incflo_level

test/test_files/abl_bndry_output/abl_bndry_output.inp

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -49,7 +49,7 @@ ABL.surface_roughness_z0 = 0.01
4949
ABL.bndry_file = "bndry_file.nc"
5050
ABL.bndry_io_mode = 0
5151
ABL.bndry_planes = ylo xlo
52-
ABL.bndry_output_start_time = 2.0
52+
ABL.bndry_output_start_time = 0.0
5353
ABL.bndry_var_names = velocity temperature
5454
ABL.bndry_output_format = netcdf
5555
#¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨#

0 commit comments

Comments
 (0)