Skip to content

Commit

Permalink
Made _nonlinear_implicit_system in TimeIntegrator non-const
Browse files Browse the repository at this point in the history
Changes:
* Made _nonlinear_implicit_system in TimeIntegrator non-const.
  In explicit time integrators, update() is called on these
  objects, so it cannot be const.
* Added doxygen for TimeIntegrator::computeTimeDerivatives()
* Fixed typo in doxygen for _du_dot_du for TimeIntegrator
* Had ActuallyExplicitEuler use _nonlinear_implicit_system
  from TimeIntegrator

Refs idaholab#11855
  • Loading branch information
joshuahansel authored and shoheiogawa committed Oct 3, 2019
1 parent 141031e commit 9a47424
Show file tree
Hide file tree
Showing 3 changed files with 28 additions and 21 deletions.
15 changes: 12 additions & 3 deletions framework/include/timeintegrators/TimeIntegrator.h
Original file line number Diff line number Diff line change
Expand Up @@ -97,6 +97,15 @@ class TimeIntegrator : public MooseObject, public Restartable
virtual void postStep() {}

virtual int order() = 0;

/**
* Computes the time derivative and the Jacobian of the time derivative
*
* This function is responsible for the following:
* - computing the time derivative; a reference is retrieved from _sys.solutionUDot().
* - computing the time derivative Jacobian, stored in _du_dot_du, which is a
* reference to _sys.duDotDu().
*/
virtual void computeTimeDerivatives() = 0;

/**
Expand Down Expand Up @@ -130,9 +139,10 @@ class TimeIntegrator : public MooseObject, public Restartable
NonlinearSystemBase & _nl;

/// Nonlinear implicit system, if applicable; otherwise, nullptr
const NonlinearImplicitSystem * _nonlinear_implicit_system;
// NonlinearImplicitSystem * _nonlinear_implicit_system;
NonlinearImplicitSystem * _nonlinear_implicit_system;

/// solution vector for \f$ {du^dot}\over{du} \f$
/// Derivative of time derivative with respect to current solution: \f$ {du^dot}\over{du} \f$
Real & _du_dot_du;
/// solution vectors
const NumericVector<Number> * const & _solution;
Expand All @@ -154,4 +164,3 @@ class TimeIntegrator : public MooseObject, public Restartable
/// Total number of linear iterations over all stages of the time step
unsigned int _n_linear_iterations;
};

32 changes: 15 additions & 17 deletions framework/src/timeintegrators/ActuallyExplicitEuler.C
Original file line number Diff line number Diff line change
Expand Up @@ -108,7 +108,7 @@ ActuallyExplicitEuler::computeTimeDerivatives()
{
if (!_sys.solutionUDot())
mooseError("ActuallyExplicitEuler: Time derivative of solution (`u_dot`) is not stored. Please "
"set uDotRequested() to true in FEProblemBase befor requesting `u_dot`.");
"set uDotRequested() to true in FEProblemBase before requesting `u_dot`.");

NumericVector<Number> & u_dot = *_sys.solutionUDot();
u_dot = *_solution;
Expand All @@ -131,29 +131,26 @@ ActuallyExplicitEuler::solve()

auto & nonlinear_system = _fe_problem.getNonlinearSystemBase();

auto & libmesh_system = dynamic_cast<NonlinearImplicitSystem &>(nonlinear_system.system());

auto & mass_matrix = *libmesh_system.matrix;
auto & mass_matrix = *_nonlinear_implicit_system->matrix;

_current_time = _fe_problem.time();

// Set time back so that we're evaluating the interior residual at the old time
_fe_problem.time() = _fe_problem.timeOld();

libmesh_system.update();
_nonlinear_implicit_system->update();

// Must compute the residual first
_explicit_residual.zero();
_fe_problem.computeResidual(*libmesh_system.current_local_solution, _explicit_residual);
_fe_problem.computeResidual(*_nonlinear_implicit_system->current_local_solution,
_explicit_residual);

// The residual is on the RHS
_explicit_residual *= -1.;

// Compute the mass matrix
_fe_problem.computeJacobianTag(*libmesh_system.current_local_solution, mass_matrix, _Ke_time_tag);

// Still testing whether leaving the old update is a good idea or not
// _explicit_euler_update = 0;
_fe_problem.computeJacobianTag(
*_nonlinear_implicit_system->current_local_solution, mass_matrix, _Ke_time_tag);

auto converged = false;

Expand Down Expand Up @@ -217,18 +214,19 @@ ActuallyExplicitEuler::solve()
mooseError("Unknown solve_type in ActuallyExplicitEuler ");
}

*libmesh_system.solution = nonlinear_system.solutionOld();
*libmesh_system.solution += _explicit_euler_update;
*_nonlinear_implicit_system->solution = nonlinear_system.solutionOld();
*_nonlinear_implicit_system->solution += _explicit_euler_update;

// Enforce contraints on the solution
DofMap & dof_map = libmesh_system.get_dof_map();
dof_map.enforce_constraints_exactly(libmesh_system, libmesh_system.solution.get());
DofMap & dof_map = _nonlinear_implicit_system->get_dof_map();
dof_map.enforce_constraints_exactly(*_nonlinear_implicit_system,
_nonlinear_implicit_system->solution.get());

libmesh_system.update();
_nonlinear_implicit_system->update();

nonlinear_system.setSolution(*libmesh_system.current_local_solution);
nonlinear_system.setSolution(*_nonlinear_implicit_system->current_local_solution);

libmesh_system.nonlinear_solver->converged = converged;
_nonlinear_implicit_system->nonlinear_solver->converged = converged;
}

void
Expand Down
2 changes: 1 addition & 1 deletion framework/src/timeintegrators/TimeIntegrator.C
Original file line number Diff line number Diff line change
Expand Up @@ -30,7 +30,7 @@ TimeIntegrator::TimeIntegrator(const InputParameters & parameters)
_fe_problem(*getCheckedPointerParam<FEProblemBase *>("_fe_problem_base")),
_sys(*getCheckedPointerParam<SystemBase *>("_sys")),
_nl(_fe_problem.getNonlinearSystemBase()),
_nonlinear_implicit_system(dynamic_cast<const NonlinearImplicitSystem *>(&_sys.system())),
_nonlinear_implicit_system(dynamic_cast<NonlinearImplicitSystem *>(&_sys.system())),
_du_dot_du(_sys.duDotDu()),
_solution(_sys.currentSolution()),
_solution_old(_sys.solutionOld()),
Expand Down

0 comments on commit 9a47424

Please sign in to comment.