From 1e7c2d3b7c2b9aa436aa07e8cac179b81c285154 Mon Sep 17 00:00:00 2001 From: Lucas Esclapez <13371051+esclapez@users.noreply.github.com> Date: Wed, 26 Jan 2022 15:25:19 -0800 Subject: [PATCH] Fix turbInflow (#47) * Fix parsing tagging input for box. * Couple of fixes for turbinflow for error introduced while re-arranging the level data. * Fix TurbInflow case to match changes in turbinflow velocity field handling. * Add a input file with refinement on loZ. --- Exec/RegTests/TurbInflow/input.3d_BoxLoZ | 78 ++++++++++++++++++++++++ Exec/RegTests/TurbInflow/pelelm_prob.H | 7 +-- Source/PeleLM.H | 2 +- Source/PeleLMBCfill.H | 30 +++++---- Source/PeleLMSetup.cpp | 1 + 5 files changed, 97 insertions(+), 21 deletions(-) create mode 100644 Exec/RegTests/TurbInflow/input.3d_BoxLoZ diff --git a/Exec/RegTests/TurbInflow/input.3d_BoxLoZ b/Exec/RegTests/TurbInflow/input.3d_BoxLoZ new file mode 100644 index 00000000..178d9899 --- /dev/null +++ b/Exec/RegTests/TurbInflow/input.3d_BoxLoZ @@ -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 diff --git a/Exec/RegTests/TurbInflow/pelelm_prob.H b/Exec/RegTests/TurbInflow/pelelm_prob.H index dd4e8c57..8aac7143 100644 --- a/Exec/RegTests/TurbInflow/pelelm_prob.H +++ b/Exec/RegTests/TurbInflow/pelelm_prob.H @@ -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; diff --git a/Source/PeleLM.H b/Source/PeleLM.H index d880a749..b73c6dbd 100644 --- a/Source/PeleLM.H +++ b/Source/PeleLM.H @@ -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 diff --git a/Source/PeleLMBCfill.H b/Source/PeleLMBCfill.H index 172e8c84..2fc3edea 100644 --- a/Source/PeleLMBCfill.H +++ b/Source/PeleLMBCfill.H @@ -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 @@ -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 @@ -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 } - } } } diff --git a/Source/PeleLMSetup.cpp b/Source/PeleLMSetup.cpp index 197be33e..c1f0a468 100644 --- a/Source/PeleLMSetup.cpp +++ b/Source/PeleLMSetup.cpp @@ -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()); }