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

more fixes for amr.subcycling_mode = None + retry #2824

Merged
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
46 changes: 23 additions & 23 deletions Exec/science/wdmerger/ci-benchmarks/wdmerger_collision_2D.out
Original file line number Diff line number Diff line change
@@ -1,29 +1,29 @@
plotfile = plt00095
time = 1.25
variables minimum value maximum value
density 8.7135211913e-05 13353524.435
xmom -4.4634547423e+14 1.4982445308e+15
ymom -1.8881098335e+15 1.9802028197e+15
density 8.6596894178e-05 14342742.997
xmom -5.3795431399e+14 1.5534999809e+15
ymom -1.8900763787e+15 1.9125216773e+15
zmom 0 0
rho_E 7.6866508899e+11 5.6480764484e+24
rho_e 7.3273369729e+11 5.6257109543e+24
Temp 100000 3970499251.5
rho_He4 8.7135211913e-17 1472.557835
rho_C12 3.4854084765e-05 5034202.6702
rho_O16 5.2281127147e-05 7781357.4533
rho_Ne20 8.7135211913e-17 643712.05721
rho_Mg24 8.7135211913e-17 1137247.8977
rho_Si28 8.7135211913e-17 4245283.1554
rho_S32 8.7135211913e-17 2178470.1687
rho_Ar36 8.7135211913e-17 497477.29422
rho_Ca40 8.7135211913e-17 382010.10154
rho_Ti44 8.7135211913e-17 1577.2314864
rho_Cr48 8.7135211913e-17 1468.2939271
rho_Fe52 8.7135211913e-17 14839.719434
rho_Ni56 8.7135211913e-17 182888.36194
phiGrav -4.614866944e+17 -2.2055818115e+16
grav_x -461232534.93 -48603.543258
grav_y -444759175.11 392332867.05
rho_E 1.0478130927e+12 1.1816191282e+25
rho_e 9.9978448556e+11 1.1805863756e+25
Temp 100000 4617810497.4
rho_He4 8.6596894178e-17 36171.490906
rho_C12 3.4638757671e-05 5329342.135
rho_O16 5.1958136506e-05 8013661.8633
rho_Ne20 8.6596894178e-17 641339.21986
rho_Mg24 8.6596894178e-17 885671.32868
rho_Si28 8.6596894178e-17 7306845.0691
rho_S32 8.6596894178e-17 4315180.4037
rho_Ar36 8.6596894178e-17 1186734.5033
rho_Ca40 8.6596894178e-17 860488.91828
rho_Ti44 8.6596894178e-17 4978.4502978
rho_Cr48 8.6596894178e-17 17110.643453
rho_Fe52 8.6596894178e-17 110923.97673
rho_Ni56 8.6596894178e-17 827475.12913
phiGrav -4.6592040725e+17 -2.205539637e+16
grav_x -466332248.65 -48558.246229
grav_y -470363988.33 399464334.41
grav_z 0 0
rho_enuc -7.502262752e+21 6.4140001008e+26
rho_enuc -4.0195280523e+22 2.0312948109e+26

124 changes: 71 additions & 53 deletions Source/driver/Castro_advance_ctu.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -157,11 +157,18 @@ Castro::retry_advance_ctu(Real dt, const advance_status& status)

if (do_retry) {

if (status.suggested_dt > 0.0_rt && status.suggested_dt < dt) {
dt_subcycle = status.suggested_dt;
int max_level_to_advance = level;

if (parent->subcyclingMode() == "None" && level == 0) {
max_level_to_advance = parent->finestLevel();
}
else {
dt_subcycle = std::min(dt, dt_subcycle) * retry_subcycle_factor;

for (int lev = level; lev <= max_level_to_advance; ++lev) {
if (status.suggested_dt > 0.0_rt && status.suggested_dt < dt) {
getLevel(lev).dt_subcycle = status.suggested_dt;
} else {
getLevel(lev).dt_subcycle = std::min(dt, getLevel(lev).dt_subcycle) * retry_subcycle_factor;
}
}

if (verbose && ParallelDescriptor::IOProcessor()) {
Expand All @@ -178,53 +185,56 @@ Castro::retry_advance_ctu(Real dt, const advance_status& status)
// be useful to us at the end of the timestep when we need
// to restore the original old data.

save_data_for_retry();
for (int lev = level; lev <= max_level_to_advance; ++lev) {
getLevel(lev).save_data_for_retry();

// Clear the contribution to the fluxes from this step.
// Clear the contribution to the fluxes from this step.

for (int dir = 0; dir < 3; ++dir) {
fluxes[dir]->setVal(0.0);
}
for (int dir = 0; dir < 3; ++dir) {
getLevel(lev).fluxes[dir]->setVal(0.0);
}

for (int dir = 0; dir < 3; ++dir) {
mass_fluxes[dir]->setVal(0.0);
}
for (int dir = 0; dir < 3; ++dir) {
getLevel(lev).mass_fluxes[dir]->setVal(0.0);
}

#if (AMREX_SPACEDIM <= 2)
if (!Geom().IsCartesian()) {
P_radial.setVal(0.0);
}
if (!Geom().IsCartesian()) {
getLevel(lev).P_radial.setVal(0.0);
}
#endif

#ifdef RADIATION
if (Radiation::rad_hydro_combined) {
for (int dir = 0; dir < AMREX_SPACEDIM; ++dir) {
rad_fluxes[dir]->setVal(0.0);
}
}
if (Radiation::rad_hydro_combined) {
for (int dir = 0; dir < AMREX_SPACEDIM; ++dir) {
getLevel(lev).rad_fluxes[dir]->setVal(0.0);
}
}
#endif

#ifdef REACTIONS
if (castro::store_burn_weights) {
burn_weights.setVal(0.0);
}
if (castro::store_burn_weights) {
getLevel(lev).burn_weights.setVal(0.0);
}
#endif

// For simplified SDC, we'll have garbage data if we
// attempt to use the lagged source terms (both reacting
// and non-reacting) from the last timestep, since that
// advance failed and we don't know if we can trust it.
// So we zero out both source term correctors.
// For simplified SDC, we'll have garbage data if we
// attempt to use the lagged source terms (both reacting
// and non-reacting) from the last timestep, since that
// advance failed and we don't know if we can trust it.
// So we zero out both source term correctors.

if (time_integration_method == SimplifiedSpectralDeferredCorrections) {
source_corrector.setVal(0.0, source_corrector.nGrow());
if (time_integration_method == SimplifiedSpectralDeferredCorrections) {
getLevel(lev).source_corrector.setVal(0.0, getLevel(lev).source_corrector.nGrow());

#ifdef SIMPLIFIED_SDC
#ifdef REACTIONS
MultiFab& SDC_react_new = get_new_data(Simplified_SDC_React_Type);
SDC_react_new.setVal(0.0, SDC_react_new.nGrow());
MultiFab& SDC_react_new = getLevel(lev).get_new_data(Simplified_SDC_React_Type);
SDC_react_new.setVal(0.0, SDC_react_new.nGrow());
#endif
#endif
}

}

}
Expand Down Expand Up @@ -311,7 +321,9 @@ Castro::subcycle_advance_ctu(const Real time, const Real dt, int amr_iteration,

if (num_subcycles_remaining > max_subcycles) {
amrex::Print() << std::endl
<< " The subcycle mechanism requested " << num_subcycles_remaining << " subcycled timesteps, which is larger than the maximum of " << max_subcycles << "." << std::endl
<< " The subcycle mechanism requested " << num_subcycles_remaining
<< " subcycled timesteps, which is larger than the maximum of "
<< max_subcycles << "." << std::endl
<< " If you would like to override this, increase the parameter castro.max_subcycles." << std::endl;
amrex::Abort("Error: too many subcycles.");
}
Expand All @@ -320,9 +332,12 @@ Castro::subcycle_advance_ctu(const Real time, const Real dt, int amr_iteration,

if (verbose && ParallelDescriptor::IOProcessor()) {
std::cout << std::endl;
std::cout << Font::Bold << FGColor::Green << " Beginning subcycle " << sub_iteration + 1 << " starting at time " << subcycle_time
std::cout << Font::Bold << FGColor::Green
<< " Beginning subcycle " << sub_iteration + 1
<< " starting at time " << subcycle_time
<< " with dt = " << dt_subcycle << ResetDisplay << std::endl;
std::cout << " Estimated number of subcycles remaining: " << num_subcycles_remaining << std::endl << std::endl;
std::cout << " Estimated number of subcycles remaining: "
<< num_subcycles_remaining << std::endl << std::endl;
}

// Swap the time levels. Only do this after the first iteration;
Expand All @@ -331,19 +346,18 @@ Castro::subcycle_advance_ctu(const Real time, const Real dt, int amr_iteration,

if (do_swap) {

swap_state_time_levels(0.0);
for (int lev = level; lev <= max_level_to_advance; ++lev) {
getLevel(lev).swap_state_time_levels(0.0);

#ifdef GRAVITY
if (do_grav) {
gravity->swapTimeLevels(level);
}
if (do_grav) {
getLevel(lev).gravity->swapTimeLevels(lev);
}
#endif
}

}
else {

} else {
do_swap = true;

}

// Set the relevant time levels.
Expand Down Expand Up @@ -462,7 +476,9 @@ Castro::subcycle_advance_ctu(const Real time, const Real dt, int amr_iteration,

// Record the number of subcycles we took for diagnostic purposes.

num_subcycles_taken = sub_iteration;
for (int lev = level; lev <= max_level_to_advance; ++lev) {
getLevel(lev).num_subcycles_taken = sub_iteration;
}

if (sub_iteration > 1) {

Expand All @@ -472,14 +488,14 @@ Castro::subcycle_advance_ctu(const Real time, const Real dt, int amr_iteration,
// we still have the last iteration's old data if we need
// it later.

for (int k = 0; k < num_state_type; k++) {

if (prev_state[k]->hasOldData()) {
state[k].replaceOldData(*prev_state[k]);
for (int lev = level; lev <= max_level_to_advance; ++lev) {
for (int k = 0; k < num_state_type; k++) {
if (getLevel(lev).prev_state[k]->hasOldData()) {
getLevel(lev).state[k].replaceOldData(*getLevel(lev).prev_state[k]);
}
getLevel(lev).state[k].setTimeLevel(time + dt, dt, 0.0);
getLevel(lev).prev_state[k]->setTimeLevel(time + dt, dt_subcycle, 0.0);
}
state[k].setTimeLevel(time + dt, dt, 0.0);
prev_state[k]->setTimeLevel(time + dt, dt_subcycle, 0.0);

}

// If we took more than one step and are going to do a reflux,
Expand All @@ -491,8 +507,10 @@ Castro::subcycle_advance_ctu(const Real time, const Real dt, int amr_iteration,
// reflux immediately following this, skip this if we're on the
// finest level and this is not the last iteration.

if (!(amr_iteration < amr_ncycle && level == parent->finestLevel())) {
keep_prev_state = true;
for (int lev = level; lev <= max_level_to_advance; ++lev) {
if (!(amr_iteration < amr_ncycle && lev == parent->finestLevel())) {
getLevel(lev).keep_prev_state = true;
}
}

}
Expand Down