Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Fix turbInflow #47

Merged
merged 5 commits into from
Jan 26, 2022
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
78 changes: 78 additions & 0 deletions Exec/RegTests/TurbInflow/input.3d_BoxLoZ
Original file line number Diff line number Diff line change
@@ -0,0 +1,78 @@
#----------------------DOMAIN DEFINITION------------------------
geometry.is_periodic = 1 1 0 # For each dir, 0: non-perio, 1: periodic
geometry.coord_sys = 0 # 0 => cart, 1 => RZ
geometry.prob_lo = 0.0 0.0 0.0 # x_lo y_lo (z_lo)
geometry.prob_hi = 0.01 0.01 0.01 # x_hi y_hi (z_hi)

# >>>>>>>>>>>>> BC FLAGS <<<<<<<<<<<<<<<<
# Interior, Inflow, Outflow, Symmetry,
# SlipWallAdiab, NoSlipWallAdiab, SlipWallIsotherm, NoSlipWallIsotherm
peleLM.lo_bc = Interior Interior Inflow
peleLM.hi_bc = Interior Interior Outflow


#-------------------------AMR CONTROL----------------------------
amr.n_cell = 32 32 32 # Level 0 number of cells in each direction
amr.v = 1 # AMR verbose
amr.max_level = 1 # maximum level number allowed
amr.ref_ratio = 2 2 2 2 # refinement ratio
amr.regrid_int = 2 # how often to regrid
amr.n_error_buf = 0 1 2 2 # number of buffer cells in error est
amr.grid_eff = 0.7 # what constitutes an efficient grid
amr.blocking_factor = 16 # block factor in grid generation (min box size)
amr.max_grid_size = 64 # max box size


#--------------------------- Problem -------------------------------
prob.T_mean = 298.0
prob.P_mean = 101325.0
prob.flowDir = 2
prob.flowMag = 5

prob.turb_file = TurbFileHIT/TurbTEST
prob.turb_scale_loc = 633.151
prob.turb_scale_vel = 0.5
prob.turb_center = 0.005 0.005
prob.turb_conv_vel = 5.
prob.turb_nplane = 32

#------------ INPUTS TO CONSTANT TRANSPORT -----------------
transport.const_viscosity = 0.0001
transport.const_bulk_viscosity = 0.0
transport.const_conductivity = 0.0
transport.const_diffusivity = 0.0001

#-------------------------PeleLM CONTROL----------------------------
peleLM.v = 1
peleLM.incompressible = 0
peleLM.rho = 1.17
peleLM.mu = 0.0
peleLM.sdc_iterMax = 1
peleLM.floor_species = 0
peleLM.num_divu_iter = 1
peleLM.num_init_iter = 3

peleLM.do_temporals = 1
peleLM.temporal_int = 2
peleLM.mass_balance = 1

amr.plot_int = 10
amr.max_step = 300
amr.dt_shrink = 0.1
amr.stop_time = 0.02
amr.cfl = 0.5
amr.derive_plot_vars = avg_pressure mag_vort mass_fractions

#--------------------REFINEMENT CONTROL------------------------
amr.refinement_indicators = loBoxZ
#amr.O2.max_level = 1
#amr.O2.value_greater = 0.4
#amr.O2.field_name = Y(O2)

amr.loBoxZ.in_box_lo = 0.0 0.0 0.0
amr.loBoxZ.in_box_hi = 0.01 0.01 0.005

fabarray.mfiter_tile_size = 1024 1024 1024
#amrex.fpe_trap_invalid = 1
#amrex.fpe_trap_zero = 1
#amrex.fpe_trap_overflow = 1
7 changes: 3 additions & 4 deletions Exec/RegTests/TurbInflow/pelelm_prob.H
Original file line number Diff line number Diff line change
Expand Up @@ -105,10 +105,9 @@ bcnormal(

auto eos = pele::physics::PhysicsType::eos();

AMREX_D_TERM(s_ext[VELX] = 0.0;,
s_ext[VELY] = 0.0;,
s_ext[VELZ] = 0.0);
s_ext[prob_parm.meanFlowDir] = prob_parm.meanFlowMag;
// When using PELE_USE_TURBINFLOW, s_ext comes in with
// turbulence data. Only add the mean flow here.
s_ext[prob_parm.meanFlowDir] += prob_parm.meanFlowMag;

massfrac[O2_ID] = 0.233;
massfrac[N2_ID] = 0.767;
Expand Down
2 changes: 1 addition & 1 deletion Source/PeleLM.H
Original file line number Diff line number Diff line change
Expand Up @@ -492,7 +492,7 @@ class PeleLM : public amrex::AmrCore {
void setInflowBoundaryVel (amrex::MultiFab &a_vel, int lev, PeleLM::TimeStamp a_time);

#ifdef PELE_USE_TURBINFLOW
void fillTurbInflow(amrex::MultiFab &a_vel, int lev, const amrex::Real a_time);
void fillTurbInflow(amrex::MultiFab &a_vel, int vel_comp, int lev, const amrex::Real a_time);
#endif

// Average down operations
Expand Down
30 changes: 14 additions & 16 deletions Source/PeleLMBCfill.H
Original file line number Diff line number Diff line change
Expand Up @@ -56,9 +56,7 @@ struct PeleLMCCFillExtDirState
#ifdef PELE_USE_TURBINFLOW
// If using TurbInflow, pass in the turbulent data. User can overwrite if needed
if ( n >= VELX && n < DENSITY ) {
for (int n = 0; n < AMREX_SPACEDIM; n++) {
s_ext[VELX+n] = vel(iv,n);
}
s_ext[VELX+n] = state(iv,VELX+n);
}
#endif
// bcnormal() is defined in pelelm_prob.H in problem directory in /Exec
Expand All @@ -81,9 +79,7 @@ struct PeleLMCCFillExtDirState
#ifdef PELE_USE_TURBINFLOW
// If using TurbInflow, pass in the turbulent data. User can overwrite if needed
if ( n >= VELX && n < DENSITY ) {
for (int n = 0; n < AMREX_SPACEDIM; n++) {
s_ext[VELX+n] = vel(iv,n);
}
s_ext[VELX+n] = state(iv,VELX+n);
}
#endif
// bcnormal() is defined in pelelm_prob.H in problem directory in /Exec
Expand Down Expand Up @@ -139,31 +135,33 @@ struct PeleLMCCFillExtDirVel
// Low
if ((bc[idir] == amrex::BCType::ext_dir) and (iv[idir] < domlo[idir])) {

#ifdef PELE_USE_TURBINFLOW
// If using TurbInflow, pass in the turbulent data. User can overwrite if needed
for (int n = 0; n < AMREX_SPACEDIM; n++) {
s_ext[VELX+n] = vel(iv,VELX+n);
}
#endif
// bcnormal() is defined in pelelm_prob.H in problem directory in /Exec
bcnormal(x, m_nAux, s_ext, idir, 1, time, geom, *lprobparm, lpmfdata);

for (int n = 0; n < AMREX_SPACEDIM; n++) {
#ifdef PELE_USE_TURBINFLOW
vel(iv,n) += s_ext[VELX+n];
#else
vel(iv,n) = s_ext[VELX+n];
#endif
}

// High
} else if ((bc[idir+AMREX_SPACEDIM] == amrex::BCType::ext_dir) and (iv[idir] > domhi[idir])) {

#ifdef PELE_USE_TURBINFLOW
// If using TurbInflow, pass in the turbulent data. User can overwrite if needed
for (int n = 0; n < AMREX_SPACEDIM; n++) {
s_ext[VELX+n] = vel(iv,VELX+n);
}
#endif
// bcnormal() is defined in pelelm_prob.H in problem directory in /Exec
bcnormal(x, m_nAux, s_ext, idir, -1, time, geom, *lprobparm, lpmfdata);

for (int n = 0; n < AMREX_SPACEDIM; n++) {
#ifdef PELE_USE_TURBINFLOW
vel(iv,n) += s_ext[VELX+n];
#else
vel(iv,n) = s_ext[VELX+n];
#endif
}

}
}
}
Expand Down
1 change: 1 addition & 0 deletions Source/PeleLMSetup.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -592,6 +592,7 @@ void PeleLM::taggingSetup()
itexists = derive_lst.canDerive(field) || isStateVariable(field);
} else if (realbox.ok()) {
errTags.push_back(AMRErrorTag(info));
itexists = true;
} else {
Abort(std::string("Unrecognized refinement indicator for " + refinement_indicators[n]).c_str());
}
Expand Down