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 restart and input_sounding issues. #1341

Merged
merged 1 commit into from
Dec 12, 2023
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
26 changes: 9 additions & 17 deletions Source/IO/Checkpoint.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -224,7 +224,7 @@ ERF::ReadCheckpointFile ()
{
amrex::Print() << "Restart from checkpoint " << restart_chkfile << "\n";

int ncomp_cons = vars_new[0][Vars::cons].nComp();


// Header
std::string File(restart_chkfile + "/Header");
Expand All @@ -238,7 +238,7 @@ ERF::ReadCheckpointFile ()

std::string line, word;

int chk_ncomp;
int chk_ncomp_cons, chk_ncomp;

// read in title line
std::getline(is, line);
Expand All @@ -251,16 +251,8 @@ ERF::ReadCheckpointFile ()
// for each variable we store

// conservative, cell-centered vars
is >> chk_ncomp;
is >> chk_ncomp_cons;
GotoNextLine(is);
#if 0
if (solverChoice.moisture_type != MoistureType::None) {
AMREX_ASSERT(chk_ncomp == ncomp_cons);
} else {
// This assumes that we are carrying RhoQ1 and RhoQ2 but not actually using them
AMREX_ASSERT(chk_ncomp == ncomp_cons-2);
}
#endif

// x-velocity on faces
is >> chk_ncomp;
Expand Down Expand Up @@ -320,18 +312,18 @@ ERF::ReadCheckpointFile ()
MakeNewLevelFromScratch (lev, t_new[lev], ba, dm);
}

// ncomp is only valid after we MakeNewLevelFromScratch (asks micro how many vars)
// NOTE: Data is written over ncomp, so check that we match the header file
int ncomp_cons = vars_new[0][Vars::cons].nComp();
AMREX_ASSERT(chk_ncomp_cons == ncomp_cons);

// read in the MultiFab data
for (int lev = 0; lev <= finest_level; ++lev)
{
MultiFab cons(grids[lev],dmap[lev],ncomp_cons,0);
VisMF::Read(cons, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "Cell"));

if (solverChoice.moisture_type != MoistureType::None) {
MultiFab::Copy(vars_new[lev][Vars::cons],cons,0,0,ncomp_cons,0);
} else {
// This assumes that we are carrying RhoQ1 and RhoQ2 but not actually using them
MultiFab::Copy(vars_new[lev][Vars::cons],cons,0,0,ncomp_cons-2,0);
}
MultiFab::Copy(vars_new[lev][Vars::cons],cons,0,0,ncomp_cons,0);

MultiFab xvel(convert(grids[lev],IntVect(1,0,0)),dmap[lev],1,0);
VisMF::Read(xvel, amrex::MultiFabFileFullPrefix(lev, restart_chkfile, "Level_", "XFace"));
Expand Down
39 changes: 16 additions & 23 deletions Source/Initialization/ERF_init_from_input_sounding.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -14,6 +14,7 @@ void
init_bx_scalars_from_input_sounding (const amrex::Box &bx,
amrex::Array4<amrex::Real> const &state,
amrex::GeometryData const &geomdata,
const bool& l_moist,
InputSoundingData const &inputSoundingData);
void
init_bx_scalars_from_input_sounding_hse (const amrex::Box &bx,
Expand All @@ -22,7 +23,9 @@ init_bx_scalars_from_input_sounding_hse (const amrex::Box &bx,
amrex::Array4<amrex::Real> const &p_hse_arr,
amrex::Array4<amrex::Real> const &pi_hse_arr,
amrex::GeometryData const &geomdata,
amrex::Real l_gravity, amrex::Real l_rdOcp,
const amrex::Real& l_gravity,
const amrex::Real& l_rdOcp,
const bool& l_moist,
InputSoundingData const &inputSoundingData);

void
Expand Down Expand Up @@ -61,6 +64,7 @@ ERF::init_from_input_sounding (int lev)

const Real l_gravity = solverChoice.gravity;
const Real l_rdOcp = solverChoice.rdOcp;
const bool l_moist = (solverChoice.moisture_type != MoistureType::None);

#ifdef _OPENMP
#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
Expand All @@ -81,14 +85,14 @@ ERF::init_from_input_sounding (int lev)
init_bx_scalars_from_input_sounding_hse(
bx, cons_arr,
r_hse_arr, p_hse_arr, pi_hse_arr,
geom[lev].data(), l_gravity, l_rdOcp, input_sounding_data);
geom[lev].data(), l_gravity, l_rdOcp, l_moist, input_sounding_data);
}
else
{
// HSE will be calculated later with call to initHSE
init_bx_scalars_from_input_sounding(
bx, cons_arr,
geom[lev].data(), input_sounding_data);
geom[lev].data(), l_moist, input_sounding_data);
}

init_bx_velocities_from_input_sounding(
Expand All @@ -111,15 +115,12 @@ void
init_bx_scalars_from_input_sounding (const amrex::Box &bx,
amrex::Array4<amrex::Real> const &state,
amrex::GeometryData const &geomdata,
const bool& l_moist,
InputSoundingData const &inputSoundingData)
{
const Real* z_inp_sound = inputSoundingData.z_inp_sound_d.dataPtr();
const Real* theta_inp_sound = inputSoundingData.theta_inp_sound_d.dataPtr();
#if defined(ERF_USE_MOISTURE)
const Real* qv_inp_sound = inputSoundingData.qv_inp_sound_d.dataPtr();
#elif defined(ERF_USE_WARM_NO_PRECIP)
const Real* qv_inp_sound = inputSoundingData.qv_inp_sound_d.dataPtr();
#endif
const int inp_sound_size = inputSoundingData.size();

// We want to set the lateral BC values, too
Expand All @@ -144,11 +145,9 @@ init_bx_scalars_from_input_sounding (const amrex::Box &bx,
state(i, j, k, RhoScalar_comp) = 0;

// total nonprecipitating water (Q1) == water vapor (Qv), i.e., there is no cloud water or cloud ice
#if defined(ERF_USE_MOISTURE)
state(i, j, k, RhoQ1_comp) = rho_0 * interpolate_1d(z_inp_sound, qv_inp_sound, z, inp_sound_size);
#elif defined(ERF_USE_WARM_NO_PRECIP)
state(i, j, k, RhoQv_comp) = rho_0 * interpolate_1d(z_inp_sound, qv_inp_sound, z, inp_sound_size);
#endif
if (l_moist)
state(i, j, k, RhoQ1_comp) = rho_0 * interpolate_1d(z_inp_sound, qv_inp_sound, z, inp_sound_size);

});
}

Expand All @@ -173,18 +172,15 @@ init_bx_scalars_from_input_sounding_hse (const amrex::Box &bx,
amrex::Array4<amrex::Real> const &p_hse_arr,
amrex::Array4<amrex::Real> const &pi_hse_arr,
amrex::GeometryData const &geomdata,
const amrex::Real l_gravity,
const amrex::Real l_rdOcp,
const amrex::Real& l_gravity,
const amrex::Real& l_rdOcp,
const bool& l_moist,
InputSoundingData const &inputSoundingData)
{
const Real* z_inp_sound = inputSoundingData.z_inp_sound_d.dataPtr();
const Real* rho_inp_sound = inputSoundingData.rho_inp_sound_d.dataPtr();
const Real* theta_inp_sound = inputSoundingData.theta_inp_sound_d.dataPtr();
#if defined(ERF_USE_MOISTURE)
const Real* qv_inp_sound = inputSoundingData.qv_inp_sound_d.dataPtr();
#elif defined(ERF_USE_WARM_NO_PRECIP)
const Real* qv_inp_sound = inputSoundingData.qv_inp_sound_d.dataPtr();
#endif
const int inp_sound_size = inputSoundingData.size();

// We want to set the lateral BC values, too
Expand Down Expand Up @@ -234,13 +230,10 @@ init_bx_scalars_from_input_sounding_hse (const amrex::Box &bx,
pi_hse_arr(i, j, k+1) = getExnergivenP(p_hse_arr(i, j, k+1), l_rdOcp);
}

#if defined(ERF_USE_MOISTURE)
// total nonprecipitating water (Q1) == water vapor (Qv), i.e., there
// is no cloud water or cloud ice
state(i, j, k, RhoQ1_comp) = rho_k * interpolate_1d(z_inp_sound, qv_inp_sound, z, inp_sound_size);
#elif defined(ERF_USE_WARM_NO_PRECIP)
state(i, j, k, RhoQv_comp) = rho_k * interpolate_1d(z_inp_sound, qv_inp_sound, z, inp_sound_size);
#endif
if (l_moist)
state(i, j, k, RhoQ1_comp) = rho_k * interpolate_1d(z_inp_sound, qv_inp_sound, z, inp_sound_size);
});
}

Expand Down
Loading