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 some small issues #1873

Merged
merged 2 commits into from
Oct 11, 2024
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
49 changes: 15 additions & 34 deletions Source/BoundaryConditions/ERF_FillPatch.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -43,20 +43,9 @@ ERF::FillPatch (int lev, Real time,
FPr_c[lev-1].FillSet(*mfs_vel[Vars::cons], time, null_bc, domain_bcs_type);
}
if (cf_set_width >= 0 && !cons_only) {
//
// This is an optimization since we won't need more than one ghost
// cell of momentum in the integrator if not using NumDiff
//
//IntVect ngu = (solverChoice.use_NumDiff) ? IntVect(1,1,1) : mfs_vel[Vars::xvel]->nGrowVect();
//IntVect ngv = (solverChoice.use_NumDiff) ? IntVect(1,1,1) : mfs_vel[Vars::yvel]->nGrowVect();
//IntVect ngw = (solverChoice.use_NumDiff) ? IntVect(1,1,0) : mfs_vel[Vars::zvel]->nGrowVect();
IntVect ngu = IntVect::TheZeroVector();
IntVect ngv = IntVect::TheZeroVector();
IntVect ngw = IntVect::TheZeroVector();

VelocityToMomentum(*mfs_vel[Vars::xvel], ngu,
*mfs_vel[Vars::yvel], ngv,
*mfs_vel[Vars::zvel], ngw,
VelocityToMomentum(*mfs_vel[Vars::xvel], IntVect{0},
*mfs_vel[Vars::yvel], IntVect{0},
*mfs_vel[Vars::zvel], IntVect{0},
*mfs_vel[Vars::cons],
*mfs_mom[IntVars::xmom],
*mfs_mom[IntVars::ymom],
Expand Down Expand Up @@ -545,16 +534,12 @@ ERF::FillIntermediatePatch (int lev, Real time,
// We always come in to this call with momenta so we need to leave with momenta!
// We need to make sure we convert back on all ghost cells/faces because this is
// how velocity from fine-fine copies (as well as physical and interpolated bcs) will be filled
if (!cons_only) {
IntVect ngu = mfs_vel[Vars::xvel]->nGrowVect();
IntVect ngv = mfs_vel[Vars::yvel]->nGrowVect();
IntVect ngw = mfs_vel[Vars::zvel]->nGrowVect();

if (!solverChoice.use_NumDiff) {
ngu = IntVect(1,1,1);
ngv = IntVect(1,1,1);
ngw = IntVect(1,1,1);
}
if (!cons_only)
{
IntVect ngu = (!solverChoice.use_NumDiff) ? IntVect(1,1,1) : mfs_vel[Vars::xvel]->nGrowVect();
IntVect ngv = (!solverChoice.use_NumDiff) ? IntVect(1,1,1) : mfs_vel[Vars::yvel]->nGrowVect();
IntVect ngw = (!solverChoice.use_NumDiff) ? IntVect(1,1,0) : mfs_vel[Vars::zvel]->nGrowVect();

VelocityToMomentum(*mfs_vel[Vars::xvel], ngu,
*mfs_vel[Vars::yvel], ngv,
*mfs_vel[Vars::zvel], ngw,
Expand Down Expand Up @@ -600,13 +585,9 @@ ERF::FillCoarsePatch (int lev, Real time)
// Convert velocity to momentum at the COARSE level
// ************************************************
//
IntVect ngu = IntVect(0,0,0);
IntVect ngv = IntVect(0,0,0);
IntVect ngw = IntVect(0,0,0);

VelocityToMomentum(vars_new[lev-1][Vars::xvel], ngu,
vars_new[lev-1][Vars::yvel], ngv,
vars_new[lev-1][Vars::zvel], ngw,
VelocityToMomentum(vars_new[lev-1][Vars::xvel], IntVect{0},
vars_new[lev-1][Vars::yvel], IntVect{0},
vars_new[lev-1][Vars::zvel], IntVect{0},
vars_new[lev-1][Vars::cons],
rU_new[lev-1],
rV_new[lev-1],
Expand Down Expand Up @@ -640,7 +621,7 @@ ERF::FillCoarsePatch (int lev, Real time)
// with InterpFromCoarseLevel which ASSUMES that all ghost cells have already been filled
// ************************************************************************************************
//
InterpFromCoarseLevel(rU_new[lev], ngu, IntVect(0,0,0), rU_new[lev-1], 0, 0, 1,
InterpFromCoarseLevel(rU_new[lev], IntVect{0}, IntVect{0}, rU_new[lev-1], 0, 0, 1,
geom[lev-1], geom[lev],
refRatio(lev-1), mapper_f, domain_bcs_type, BCVars::xvel_bc);

Expand All @@ -650,15 +631,15 @@ ERF::FillCoarsePatch (int lev, Real time)
// with InterpFromCoarseLevel which ASSUMES that all ghost cells have already been filled
// ************************************************************************************************
//
InterpFromCoarseLevel(rV_new[lev], ngv, IntVect(0,0,0), rV_new[lev-1], 0, 0, 1,
InterpFromCoarseLevel(rV_new[lev], IntVect{0}, IntVect{0}, rV_new[lev-1], 0, 0, 1,
geom[lev-1], geom[lev],
refRatio(lev-1), mapper_f, domain_bcs_type, BCVars::yvel_bc);

//************************************************************************************************
// Interpolate z-momentum from coarse to fine level
// with InterpFromCoarseLevel which ASSUMES that all ghost cells have already been filled
// ************************************************************************************************
InterpFromCoarseLevel(rW_new[lev], ngw, IntVect(0,0,0), rW_new[lev-1], 0, 0, 1,
InterpFromCoarseLevel(rW_new[lev], IntVect{0}, IntVect{0}, rW_new[lev-1], 0, 0, 1,
geom[lev-1], geom[lev],
refRatio(lev-1), mapper_f, domain_bcs_type, BCVars::zvel_bc);
//
Expand Down
4 changes: 2 additions & 2 deletions Source/Diffusion/ERF_ComputeStrain_N.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -208,8 +208,8 @@ ComputeStrain_N (Box bxcc, Box tbxxy, Box tbxxz, Box tbxyz, Box domain,
//***********************************************************************************
// Cell centered strains
ParallelFor(bxcc, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
tau11(i,j,k) = (u(i+1, j , k )/mf_u(i+1,j,0) - u(i, j, k)/mf_u(i,j,0))*dxInv[0]*mf_u(i,j,0)*mf_u(i,j,0);
tau22(i,j,k) = (v(i , j+1, k )/mf_v(i,j+1,0) - v(i, j, k)/mf_v(i,j,0))*dxInv[1]*mf_v(i,j,0)*mf_v(i,j,0);
tau11(i,j,k) = (u(i+1, j , k )/mf_u(i+1,j,0) - u(i, j, k)/mf_u(i,j,0))*dxInv[0]*mf_m(i,j,0)*mf_m(i,j,0);
tau22(i,j,k) = (v(i , j+1, k )/mf_v(i,j+1,0) - v(i, j, k)/mf_v(i,j,0))*dxInv[1]*mf_m(i,j,0)*mf_m(i,j,0);
tau33(i,j,k) = (w(i , j , k+1) - w(i, j, k))*dxInv[2];
});

Expand Down
7 changes: 4 additions & 3 deletions Source/Diffusion/ERF_ComputeStrain_T.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -25,6 +25,7 @@ using namespace amrex;
* @param[in] z_nd nodal array of physical z heights
* @param[in] bc_ptr container with boundary condition types
* @param[in] dxInv inverse cell size array
* @param[in] mf_m map factor at cell center
* @param[in] mf_u map factor at x-face
* @param[in] mf_v map factor at y-face
*/
Expand All @@ -38,7 +39,7 @@ ComputeStrain_T (Box bxcc, Box tbxxy, Box tbxxz, Box tbxyz, Box domain,
const Array4<const Real>& z_nd,
const Array4<const Real>& detJ,
const BCRec* bc_ptr, const GpuArray<Real, AMREX_SPACEDIM>& dxInv,
const Array4<const Real>& /*mf_m*/,
const Array4<const Real>& mf_m,
const Array4<const Real>& mf_u,
const Array4<const Real>& mf_v)
{
Expand Down Expand Up @@ -452,9 +453,9 @@ ComputeStrain_T (Box bxcc, Box tbxxy, Box tbxxz, Box tbxyz, Box domain,
met_h_zeta = detJ(i,j,k);

tau11(i,j,k) = ( (u(i+1, j, k)/mf_u(i+1,j,0) - u(i, j, k)/mf_u(i,j,0))*dxInv[0]
- (met_h_xi/met_h_zeta)*GradUz ) * mf_u(i,j,0)*mf_u(i,j,0);
- (met_h_xi/met_h_zeta)*GradUz ) * mf_m(i,j,0)*mf_m(i,j,0);
tau22(i,j,k) = ( (v(i, j+1, k)/mf_v(i,j+1,0) - v(i, j, k)/mf_v(i,j,0))*dxInv[1]
- (met_h_eta/met_h_zeta)*GradVz ) * mf_v(i,j,0)*mf_v(i,j,0);
- (met_h_eta/met_h_zeta)*GradVz ) * mf_m(i,j,0)*mf_m(i,j,0);
tau33(i,j,k) = (w(i, j, k+1) - w(i, j, k))*dxInv[2]/met_h_zeta;
});

Expand Down
3 changes: 3 additions & 0 deletions Source/ERF_make_new_arrays.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -147,6 +147,9 @@ ERF::init_stuff (int lev, const BoxArray& ba, const DistributionMapping& dm,
rW_new[lev].define(convert(ba, IntVect(0,0,1)), dm, 1, ngrow_vels);

// We do this here just so they won't be undefined in the initial FillPatch
rU_old[lev].setVal(1.2e21);
rV_old[lev].setVal(3.4e22);
rW_old[lev].setVal(5.6e23);
rU_new[lev].setVal(1.2e21);
rV_new[lev].setVal(3.4e22);
rW_new[lev].setVal(5.6e23);
Expand Down
2 changes: 1 addition & 1 deletion Source/ERF_make_new_level.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -26,7 +26,7 @@ void ERF::MakeNewLevelFromScratch (int lev, Real time, const BoxArray& ba_in,
BoxArray ba;
DistributionMapping dm;
Box domain(Geom(0).Domain());
if (lev == 0 &&
if (lev == 0 && restart_chkfile.empty() &&
(max_grid_size[0][0] >= domain.length(0)) &&
(max_grid_size[0][1] >= domain.length(1)) &&
ba_in.size() != ParallelDescriptor::NProcs())
Expand Down
7 changes: 4 additions & 3 deletions Source/TimeIntegration/ERF_advance_dycore.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -252,9 +252,10 @@ void ERF::advance_dycore(int level,
// This is an optimization since we won't need more than one ghost
// cell of momentum in the integrator if not using NumDiff
//
IntVect ngu = (solverChoice.use_NumDiff) ? IntVect(1,1,1) : xvel_old.nGrowVect();
IntVect ngv = (solverChoice.use_NumDiff) ? IntVect(1,1,1) : yvel_old.nGrowVect();
IntVect ngw = (solverChoice.use_NumDiff) ? IntVect(1,1,0) : zvel_old.nGrowVect();
IntVect ngu = (!solverChoice.use_NumDiff) ? IntVect(1,1,1) : xvel_old.nGrowVect();
IntVect ngv = (!solverChoice.use_NumDiff) ? IntVect(1,1,1) : yvel_old.nGrowVect();
IntVect ngw = (!solverChoice.use_NumDiff) ? IntVect(1,1,0) : zvel_old.nGrowVect();

VelocityToMomentum(xvel_old, ngu, yvel_old, ngv, zvel_old, ngw, density,
state_old[IntVars::xmom],
state_old[IntVars::ymom],
Expand Down
20 changes: 7 additions & 13 deletions Source/Utils/ERF_VelocityToMomentum.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -54,19 +54,9 @@ void VelocityToMomentum (const MultiFab& xvel_in,
tby = mfi.tilebox(IntVect(0,1,0),yvel_ngrow);
tbz = mfi.tilebox(IntVect(0,0,1),zvel_ngrow);

#if 0
if (l_use_ndiff) {
tbx = mfi.tilebox(IntVect(1,0,0),xvel_ngrow);
tby = mfi.tilebox(IntVect(0,1,0),yvel_ngrow);
tbz = mfi.tilebox(IntVect(0,0,1),zvel_ngrow);
} else {
tbx = mfi.tilebox(IntVect(1,0,0),IntVect(1,1,1));
if (tbx.smallEnd(2) < 0) tbx.setSmall(2,0);
tby = mfi.tilebox(IntVect(0,1,0),IntVect(1,1,1));
if (tby.smallEnd(2) < 0) tby.setSmall(2,0);
tbz = mfi.tilebox(IntVect(0,0,1),IntVect(1,1,0));
}
#endif
// Don't actually try to fill w above or below the domain
if (tbz.smallEnd(2) < domain.smallEnd(2)) tbz.setSmall(2,domain.smallEnd(2));
if (tbz.bigEnd(2) > domain.bigEnd(2)+1) tbz.setBig(2,domain.bigEnd(2)+1);

// Conserved/state variables on cell centers -- we use this for density
const Array4<const Real>& dens_arr = density.array(mfi);
Expand All @@ -81,6 +71,8 @@ void VelocityToMomentum (const MultiFab& xvel_in,
const Array4<Real const>& vely = yvel_in.const_array(mfi);
const Array4<Real const>& velz = zvel_in.const_array(mfi);

// ********************************************************************************************

ParallelFor(tbx, tby, tbz,
[=] AMREX_GPU_DEVICE (int i, int j, int k) {
momx(i,j,k) = velx(i,j,k) * 0.5 * (dens_arr(i,j,k,Rho_comp) + dens_arr(i-1,j,k,Rho_comp));
Expand All @@ -92,6 +84,8 @@ void VelocityToMomentum (const MultiFab& xvel_in,
momz(i,j,k) = velz(i,j,k) * 0.5 * (dens_arr(i,j,k,Rho_comp) + dens_arr(i,j,k-1,Rho_comp));
});

// ********************************************************************************************

if ( (bx.smallEnd(0) == domain.smallEnd(0)) &&
(bc_ptr_h[BCVars::cons_bc].lo(0) == ERFBCType::ext_dir) ) {
ParallelFor(makeSlab(tbx,0,domain.smallEnd(0)), [=] AMREX_GPU_DEVICE (int i, int j, int k) {
Expand Down
Loading