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

Add and test skip_times option for Predictor #6507

Merged
merged 2 commits into from Mar 9, 2016
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.
Jump to
Jump to file
Failed to load files.
Diff view
Diff view
3 changes: 3 additions & 0 deletions framework/include/postprocessors/Residual.h
Expand Up @@ -35,6 +35,9 @@ class Residual : public GeneralPostprocessor
* This will return the final nonlinear residual.
*/
virtual Real getValue();

protected:
MooseEnum _residual_type;
};

#endif // RESIDUAL_H
3 changes: 2 additions & 1 deletion framework/include/predictors/AdamsPredictor.h
Expand Up @@ -40,9 +40,10 @@ class AdamsPredictor : public Predictor
virtual ~AdamsPredictor();

virtual int order() { return _order; }
virtual void timestepSetup();
virtual bool shouldApply();
virtual void apply(NumericVector<Number> & sln);
virtual NumericVector<Number> & solutionPredictor() { return _solution_predictor; }
virtual void historyControl();

protected:
int _order;
Expand Down
5 changes: 5 additions & 0 deletions framework/include/predictors/Predictor.h
Expand Up @@ -44,6 +44,8 @@ class Predictor :
virtual ~Predictor();

virtual int order() { return 0; }
virtual void timestepSetup();
virtual bool shouldApply();
virtual void apply(NumericVector<Number> & sln) = 0;

virtual NumericVector<Number> & solutionPredictor() { return _solution_predictor; }
Expand All @@ -62,6 +64,9 @@ class Predictor :

/// Amount by which to scale the predicted value. Must be in [0,1].
Real _scale;

/// Times for which the predictor should not be applied
std::vector<Real> _skip_times;
};

#endif /* PREDICTOR_H */
4 changes: 1 addition & 3 deletions framework/include/predictors/SimplePredictor.h
Expand Up @@ -50,10 +50,8 @@ class SimplePredictor : public Predictor
SimplePredictor(const InputParameters & parameters);
virtual ~SimplePredictor();

virtual bool shouldApply();
virtual void apply(NumericVector<Number> & sln);

protected:
Real _scale;
};

#endif /* SIMPLEPREDICTOR_H */
4 changes: 3 additions & 1 deletion framework/src/base/NonlinearSystem.C
Expand Up @@ -784,13 +784,15 @@ void
NonlinearSystem::onTimestepBegin()
{
_time_integrator->preSolve();
if (_predictor.get())
_predictor->timestepSetup();
}

void
NonlinearSystem::setInitialSolution()
{
NumericVector<Number> & initial_solution(solution());
if (_predictor.get())
if (_predictor.get() && _predictor->shouldApply())
{
_predictor->apply(initial_solution);
_fe_problem.predictorCleanup(initial_solution);
Expand Down
24 changes: 22 additions & 2 deletions framework/src/postprocessors/Residual.C
Expand Up @@ -16,21 +16,41 @@

#include "FEProblem.h"
#include "SubProblem.h"
#include "NonlinearSystem.h"

template<>
InputParameters validParams<Residual>()
{
MooseEnum residual_types("FINAL INITIAL_BEFORE_PRESET INITIAL_AFTER_PRESET", "FINAL");

InputParameters params = validParams<GeneralPostprocessor>();
params.addParam<MooseEnum>("residual_type", residual_types, "Type of residual to be reported. Choices are: "+residual_types.getRawNames());
return params;
}

Residual::Residual(const InputParameters & parameters) :
GeneralPostprocessor(parameters)
GeneralPostprocessor(parameters),
_residual_type(getParam<MooseEnum>("residual_type"))
{}

Real
Residual::getValue()
{
return _subproblem.finalNonlinearResidual();
Real residual = 0.0;
if (_residual_type == "FINAL")
residual = _subproblem.finalNonlinearResidual();
else
{
FEProblem * fe_problem = dynamic_cast<FEProblem *> (&_subproblem);
if (!fe_problem)
mooseError("Dynamic cast to FEProblem failed in Residual Postprocessor");
if (_residual_type == "INITIAL_BEFORE_PRESET")
residual = fe_problem->getNonlinearSystem()._initial_residual_before_preset_bcs;
else if (_residual_type == "INITIAL_AFTER_PRESET")
residual = fe_problem->getNonlinearSystem()._initial_residual_after_preset_bcs;
else
mooseError("Invalid residual_type option in Residual Postprocessor: "<<_residual_type);
}
return residual;
}

21 changes: 13 additions & 8 deletions framework/src/predictors/AdamsPredictor.C
Expand Up @@ -46,7 +46,7 @@ AdamsPredictor::~AdamsPredictor()
}

void
AdamsPredictor::historyControl()
AdamsPredictor::timestepSetup()
{
// if the time step number hasn't changed then do nothing
if (_t_step == _t_step_old)
Expand All @@ -66,20 +66,25 @@ AdamsPredictor::historyControl()
_dtstorage = _dt_old;
}

void
AdamsPredictor::apply(NumericVector<Number> & sln)
bool
AdamsPredictor::shouldApply()
{
// At the moment, I don't believe there is a function in Predictor
// that gets called on time step begin.
// That means that history control must go here.
historyControl();
if (!Predictor::shouldApply())
return false;

// AB2 can only be applied if there are enough old solutions
// AB1 could potentially be used for the time step prior?
// It would be possible to do VSVO Adams, Kevin has the info
// Doing so requires a time stack of some sort....
if (_dt == 0 || _dt_old == 0 || _dt_older == 0 || _t_step < 2)
return;
return false;
else
return true;
}

void
AdamsPredictor::apply(NumericVector<Number> & sln)
{
// localize current solution to working vec
sln.localize(_solution_predictor);
// NumericVector<Number> & vector1 = _tmp_previous_solution;
Expand Down
25 changes: 24 additions & 1 deletion framework/src/predictors/Predictor.C
Expand Up @@ -25,6 +25,7 @@ InputParameters validParams<Predictor>()
{
InputParameters params = validParams<MooseObject>();
params.addRequiredParam<Real>("scale", "The scale factor for the predictor (can range from 0 to 1)");
params.addParam<std::vector<Real> >("skip_times", "Time steps for which the predictor should not be applied");

params.registerBase("Predictor");

Expand All @@ -44,7 +45,8 @@ Predictor::Predictor(const InputParameters & parameters) :
_solution_old(_nl.solutionOld()),
_solution_older(_nl.solutionOlder()),
_solution_predictor(_nl.addVector("predictor", true, GHOSTED)),
_scale(getParam<Real>("scale"))
_scale(getParam<Real>("scale")),
_skip_times(getParam<std::vector<Real> >("skip_times"))
{
if (_scale < 0.0 || _scale > 1.0)
mooseError("Input value for scale = " << _scale << " is outside of permissible range (0 to 1)");
Expand All @@ -54,3 +56,24 @@ Predictor::~Predictor()
{
}

void
Predictor::timestepSetup()
{
}

bool
Predictor::shouldApply()
{
bool should_apply = true;

const Real & current_time = _fe_problem.time();
for (unsigned int i=0; i<_skip_times.size(); ++i)
{
if (MooseUtils::absoluteFuzzyEqual(current_time, _skip_times[i]))
{
should_apply = false;
break;
}
}
return should_apply;
}
42 changes: 26 additions & 16 deletions framework/src/predictors/SimplePredictor.C
Expand Up @@ -19,39 +19,49 @@ template<>
InputParameters validParams<SimplePredictor>()
{
InputParameters params = validParams<Predictor>();
params.addRequiredParam<Real>("scale", "The scale factor for the predictor (can range from 0 to 1)");

return params;
}

SimplePredictor::SimplePredictor(const InputParameters & parameters) :
Predictor(parameters),
_scale(getParam<Real>("scale"))
Predictor(parameters)
{
}

SimplePredictor::~SimplePredictor()
{
}

bool
SimplePredictor::shouldApply()
{
bool should_apply = true;
should_apply = Predictor::shouldApply();

if (_t_step < 2 || _dt_old <= 0)
should_apply = false;

if (!should_apply)
_console << " Skipping predictor this step" << std::endl;

return should_apply;
}

void
SimplePredictor::apply(NumericVector<Number> & sln)
{
if (_dt_old > 0)
{
// Save the original stream flags
std::ios_base::fmtflags out_flags = Moose::out.flags();
// Save the original stream flags
std::ios_base::fmtflags out_flags = Moose::out.flags();

_console << " Applying predictor with scale factor = " << std::fixed << std::setprecision(2) << _scale << std::endl;
_console << " Applying predictor with scale factor = " << std::fixed << std::setprecision(2) << _scale << std::endl;

// Restore the flags
Moose::out.flags(out_flags);
// Restore the flags
Moose::out.flags(out_flags);

Real dt_adjusted_scale_factor = _scale * _dt / _dt_old;
if (dt_adjusted_scale_factor != 0.0)
{
sln *= (1.0 + dt_adjusted_scale_factor);
sln.add(-dt_adjusted_scale_factor, _solution_older);
}
Real dt_adjusted_scale_factor = _scale * _dt / _dt_old;
if (dt_adjusted_scale_factor != 0.0)
{
sln *= (1.0 + dt_adjusted_scale_factor);
sln.add(-dt_adjusted_scale_factor, _solution_older);
}
}
Binary file removed modules/solid_mechanics/tests/predictor/gold/out.e
Binary file not shown.
@@ -0,0 +1,5 @@
time,initial_residual
0,0
0.5,111278.14081815
1,104770.5728571

@@ -0,0 +1,5 @@
time,initial_residual
0,0
0.5,111278.14081815
1,3.9203621909661e-10