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 "ImplicitEuler.h" 11 : #include "NonlinearSystem.h" 12 : 13 : registerMooseObject("MooseApp", ImplicitEuler); 14 : 15 : InputParameters 16 95558 : ImplicitEuler::validParams() 17 : { 18 95558 : InputParameters params = TimeIntegrator::validParams(); 19 95558 : params.addClassDescription("Time integration using the implicit Euler method."); 20 95558 : return params; 21 0 : } 22 : 23 54206 : ImplicitEuler::ImplicitEuler(const InputParameters & parameters) : TimeIntegrator(parameters) {} 24 : 25 100088 : ImplicitEuler::~ImplicitEuler() {} 26 : 27 : void 28 4116509 : ImplicitEuler::computeTimeDerivatives() 29 : { 30 4116509 : if (!_sys.solutionUDot()) 31 0 : mooseError("ImplicitEuler: Time derivative of solution (`u_dot`) is not stored. Please set " 32 : "uDotRequested() to true in FEProblemBase befor requesting `u_dot`."); 33 : 34 4116509 : NumericVector<Number> & u_dot = *_sys.solutionUDot(); 35 4116509 : if (!_var_restriction) 36 : { 37 4099727 : u_dot = *_solution; 38 4099727 : computeTimeDerivativeHelper(u_dot, _solution_old); 39 : } 40 : else 41 : { 42 16782 : auto u_dot_sub = u_dot.get_subvector(_local_indices); 43 16782 : _solution->create_subvector(*_solution_sub, _local_indices, false); 44 16782 : _solution_old.create_subvector(*_solution_old_sub, _local_indices, false); 45 16782 : *u_dot_sub = *_solution_sub; 46 16782 : computeTimeDerivativeHelper(*u_dot_sub, *_solution_old_sub); 47 16782 : u_dot.restore_subvector(std::move(u_dot_sub), _local_indices); 48 : // Scatter info needed for ghosts 49 16782 : u_dot.close(); 50 16782 : } 51 : 52 4116509 : computeDuDotDu(); 53 4116509 : } 54 : 55 : void 56 24894910 : ImplicitEuler::computeADTimeDerivatives(ADReal & ad_u_dot, 57 : const dof_id_type & dof, 58 : ADReal & /*ad_u_dotdot*/) const 59 : { 60 24894910 : computeTimeDerivativeHelper(ad_u_dot, _solution_old(dof)); 61 24894910 : } 62 : 63 : void 64 2326125 : ImplicitEuler::postResidual(NumericVector<Number> & residual) 65 : { 66 2326125 : if (!_var_restriction) 67 : { 68 2309287 : residual += *_Re_time; 69 2309287 : residual += *_Re_non_time; 70 2309287 : residual.close(); 71 : } 72 : else 73 : { 74 16838 : auto residual_sub = residual.get_subvector(_local_indices); 75 16838 : auto re_time_sub = _Re_time->get_subvector(_local_indices); 76 16838 : auto re_non_time_sub = _Re_non_time->get_subvector(_local_indices); 77 16838 : *residual_sub += *re_time_sub; 78 16838 : *residual_sub += *re_non_time_sub; 79 16838 : residual.restore_subvector(std::move(residual_sub), _local_indices); 80 16838 : _Re_time->restore_subvector(std::move(re_time_sub), _local_indices); 81 16838 : _Re_non_time->restore_subvector(std::move(re_non_time_sub), _local_indices); 82 16838 : } 83 2326125 : } 84 : 85 : Real 86 10400 : ImplicitEuler::timeDerivativeRHSContribution(dof_id_type dof_id, 87 : const std::vector<Real> & factors) const 88 : { 89 : mooseAssert(factors.size() == numStatesRequired(), 90 : "Either too many or too few states are given!"); 91 10400 : return factors[0] * _solution_old(dof_id) / _dt; 92 : } 93 : 94 : Real 95 10400 : ImplicitEuler::timeDerivativeMatrixContribution(const Real factor) const 96 : { 97 10400 : return factor / _dt; 98 : }