Introduce parameter for time step in days when single precision should be used

in the linear solvers. Also, a parameter was introduced to toggle the use of AMG.
This commit is contained in:
Robert Kloefkorn
2016-11-21 17:18:24 +01:00
parent 720d341c76
commit feea8c1753
5 changed files with 16 additions and 5 deletions

View File

@@ -255,8 +255,8 @@ typedef Eigen::Array<double,
const bool converged = asImpl().getConvergence(timer, iteration);
const bool must_solve = (iteration < nonlinear_solver.minIter()) || (!converged);
if (must_solve) {
// enable single precision for solvers when dt is smaller then 20 days
residual_.singlePrecision = (unit::convert::to(dt, unit::day) < 20.) ;
// enable single precision for solvers when dt is smaller then maximal time step for single precision
residual_.singlePrecision = ( dt < param_.maxSinglePrecisionTimeStep_ );
// Compute the nonlinear update.
V dx = asImpl().solveJacobianSystem();
@@ -2157,7 +2157,7 @@ typedef Eigen::Array<double,
// updatePhaseCondFromPrimarVariable() logic requires active gas and oil phase.
phaseCondition_.assign(nc, PhasePresence());
return;
}
}
for (int c = 0; c < nc; ++c) {
phaseCondition_[c] = PhasePresence(); // No free phases.
phaseCondition_[c].setFreeWater(); // Not necessary for property calculation usage.

View File

@@ -19,6 +19,8 @@
#include <opm/autodiff/BlackoilModelParameters.hpp>
#include <opm/core/utility/parameters/ParameterGroup.hpp>
#include <opm/parser/eclipse/Units/Units.hpp>
namespace Opm
{
@@ -47,6 +49,8 @@ namespace Opm
tolerance_cnv_ = param.getDefault("tolerance_cnv", tolerance_cnv_);
tolerance_wells_ = param.getDefault("tolerance_wells", tolerance_wells_ );
tolerance_well_control_ = param.getDefault("tolerance_well_control", tolerance_well_control_);
maxSinglePrecisionTimeStep_ = unit::convert::from(
param.getDefault("max_single_precision_days", unit::convert::to( maxSinglePrecisionTimeStep_, unit::day) ), unit::day );
solve_welleq_initially_ = param.getDefault("solve_welleq_initially",solve_welleq_initially_);
update_equations_scaling_ = param.getDefault("update_equations_scaling", update_equations_scaling_);
compute_well_potentials_ = param.getDefault("compute_well_potentials", compute_well_potentials_);
@@ -68,6 +72,7 @@ namespace Opm
tolerance_cnv_ = 1.0e-2;
tolerance_wells_ = 1.0e-3;
tolerance_well_control_ = 1.0e-7;
maxSinglePrecisionTimeStep_ = unit::convert::from( 20.0, unit::day );
solve_welleq_initially_ = true;
update_equations_scaling_ = false;
compute_well_potentials_ = false;

View File

@@ -48,6 +48,10 @@ namespace Opm
// TODO: it might need to distinguish between rate control and pressure control later
double tolerance_well_control_;
/// Tolerance for time step in seconds where single precision can be used
/// for solving for the Jacobian
double maxSinglePrecisionTimeStep_;
/// Solve well equation initially
bool solve_welleq_initially_;

View File

@@ -206,8 +206,7 @@ namespace Opm
parallelInformation_arg.copyOwnerToAll(istlb, istlb);
#if ! HAVE_UMFPACK
const bool useAmg = false ;
if( useAmg )
if( parameters_.linear_solver_use_amg_ )
{
typedef ISTLUtility::CPRSelector< Matrix, Vector, Vector, POrComm> CPRSelectorType;
typedef typename CPRSelectorType::AMG AMG;

View File

@@ -43,6 +43,7 @@ namespace Opm
bool newton_use_gmres_;
bool require_full_sparsity_pattern_;
bool ignoreConvergenceFailure_;
bool linear_solver_use_amg_;
NewtonIterationBlackoilInterleavedParameters() { reset(); }
// read values from parameter class
@@ -59,6 +60,7 @@ namespace Opm
linear_solver_verbosity_ = param.getDefault("linear_solver_verbosity", linear_solver_verbosity_);
require_full_sparsity_pattern_ = param.getDefault("require_full_sparsity_pattern", require_full_sparsity_pattern_);
ignoreConvergenceFailure_ = param.getDefault("linear_solver_ignoreconvergencefailure", ignoreConvergenceFailure_);
linear_solver_use_amg_ = param.getDefault("linear_solver_use_amg", linear_solver_use_amg_ );
}
// set default values
@@ -71,6 +73,7 @@ namespace Opm
linear_solver_verbosity_ = 0;
require_full_sparsity_pattern_ = false;
ignoreConvergenceFailure_ = false;
linear_solver_use_amg_ = false;
}
};