Line data Source code
1 : //* This file is part of the MOOSE framework 2 : //* https://mooseframework.inl.gov 3 : //* 4 : //* All rights reserved, see COPYRIGHT for full restrictions 5 : //* https://github.com/idaholab/moose/blob/master/COPYRIGHT 6 : //* 7 : //* Licensed under LGPL 2.1, please see LICENSE for details 8 : //* https://www.gnu.org/licenses/lgpl-2.1.html 9 : 10 : #include "AdamsPredictor.h" 11 : #include "NonlinearSystem.h" 12 : 13 : #include "libmesh/numeric_vector.h" 14 : 15 : registerMooseObject("MooseApp", AdamsPredictor); 16 : 17 : InputParameters 18 14313 : AdamsPredictor::validParams() 19 : { 20 14313 : InputParameters params = Predictor::validParams(); 21 14313 : params.addClassDescription( 22 : "Implements an explicit Adams predictor based on two old solution vectors."); 23 14313 : params.addParam<int>("order", 2, "The maximum reachable order of the Adams-Bashforth Predictor"); 24 14313 : return params; 25 0 : } 26 : 27 24 : AdamsPredictor::AdamsPredictor(const InputParameters & parameters) 28 : : Predictor(parameters), 29 24 : _order(getParam<int>("order")), 30 24 : _current_old_solution(_nl.addVector("AB2_current_old_solution", true, libMesh::GHOSTED)), 31 24 : _older_solution(_nl.addVector("AB2_older_solution", true, libMesh::GHOSTED)), 32 24 : _oldest_solution(_nl.addVector("AB2_rejected_solution", true, libMesh::GHOSTED)), 33 24 : _tmp_previous_solution(_nl.addVector("tmp_previous_solution", true, libMesh::GHOSTED)), 34 24 : _tmp_residual_old(_nl.addVector("tmp_residual_old", true, libMesh::GHOSTED)), 35 24 : _tmp_third_vector(_nl.addVector("tmp_third_vector", true, libMesh::GHOSTED)), 36 24 : _dt_older(declareRestartableData<Real>("dt_older", 0)), 37 48 : _dtstorage(declareRestartableData<Real>("dtstorage", 0)) 38 : { 39 24 : } 40 : 41 : void 42 88 : AdamsPredictor::timestepSetup() 43 : { 44 88 : Predictor::timestepSetup(); 45 : 46 : // if the time step number hasn't changed (=repeated timestep) then do nothing 47 88 : if (_is_repeated_timestep) 48 0 : return; 49 : 50 : // Otherwise move back the previous old solution and copy the current old solution, 51 : // This will probably not work with DT2, but I don't need to get it to work with dt2. 52 88 : _older_solution.localize(_oldest_solution); 53 : // Set older solution to hold the previous old solution 54 88 : _current_old_solution.localize(_older_solution); 55 : // Set current old solution to hold what it says. 56 88 : (_nl.solutionOld()).localize(_current_old_solution); 57 : // Same thing for dt 58 88 : _dt_older = _dtstorage; 59 88 : _dtstorage = _dt_old; 60 : } 61 : 62 : bool 63 88 : AdamsPredictor::shouldApply() 64 : { 65 88 : if (!Predictor::shouldApply()) 66 0 : return false; 67 : 68 : // AB2 can only be applied if there are enough old solutions 69 : // AB1 could potentially be used for the time step prior? 70 : // It would be possible to do VSVO Adams, Kevin has the info 71 : // Doing so requires a time stack of some sort.... 72 88 : if (_dt == 0 || _dt_old == 0 || _dt_older == 0 || _t_step < 2) 73 22 : return false; 74 : else 75 66 : return true; 76 : } 77 : 78 : void 79 66 : AdamsPredictor::apply(NumericVector<Number> & sln) 80 : { 81 : // localize current solution to working vec 82 66 : sln.localize(_solution_predictor); 83 : // NumericVector<Number> & vector1 = _tmp_previous_solution; 84 66 : NumericVector<Number> & vector2 = _tmp_residual_old; 85 66 : NumericVector<Number> & vector3 = _tmp_third_vector; 86 : 87 66 : Real commonpart = _dt / _dt_old; 88 66 : Real firstpart = (1 + .5 * commonpart); 89 66 : Real secondpart = .5 * _dt / _dt_older; 90 : 91 66 : _older_solution.localize(vector2); 92 66 : _oldest_solution.localize(vector3); 93 : 94 66 : _solution_predictor *= 1 + commonpart * firstpart; 95 66 : vector2 *= -1. * commonpart * (firstpart + secondpart); 96 66 : vector3 *= commonpart * secondpart; 97 : 98 66 : _solution_predictor += vector2; 99 66 : _solution_predictor += vector3; 100 : 101 66 : _solution_predictor.localize(sln); 102 66 : }