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 "TimeIntegrator.h" 11 : #include "FEProblem.h" 12 : #include "NonlinearSystemBase.h" 13 : 14 : #include "libmesh/nonlinear_implicit_system.h" 15 : #include "libmesh/nonlinear_solver.h" 16 : #include "libmesh/dof_map.h" 17 : 18 : using namespace libMesh; 19 : 20 : InputParameters 21 328097 : TimeIntegrator::validParams() 22 : { 23 328097 : InputParameters params = MooseObject::validParams(); 24 328097 : params.addParam<std::vector<VariableName>>( 25 : "variables", {}, "A subset of the variables that this time integrator should be applied to"); 26 328097 : params.registerBase("TimeIntegrator"); 27 328097 : return params; 28 0 : } 29 : 30 57042 : TimeIntegrator::TimeIntegrator(const InputParameters & parameters) 31 : : MooseObject(parameters), 32 : Restartable(this, "TimeIntegrators"), 33 114084 : NonlinearTimeIntegratorInterface(*getCheckedPointerParam<FEProblemBase *>("_fe_problem_base"), 34 114084 : *getCheckedPointerParam<SystemBase *>("_sys")), 35 114084 : LinearTimeIntegratorInterface(*getCheckedPointerParam<SystemBase *>("_sys")), 36 57042 : _fe_problem(*getCheckedPointerParam<FEProblemBase *>("_fe_problem_base")), 37 57042 : _sys(*getCheckedPointerParam<SystemBase *>("_sys")), 38 57042 : _du_dot_du(_sys.duDotDus()), 39 57042 : _solution(_sys.currentSolution()), 40 57042 : _solution_old(_sys.solutionState(1)), 41 57042 : _solution_sub(declareRestartableDataWithContext<std::unique_ptr<NumericVector<Number>>>( 42 57042 : "solution_sub", &const_cast<libMesh::Parallel::Communicator &>(this->comm()))), 43 57042 : _solution_old_sub(declareRestartableDataWithContext<std::unique_ptr<NumericVector<Number>>>( 44 57042 : "solution_old_sub", &const_cast<libMesh::Parallel::Communicator &>(this->comm()))), 45 57042 : _t_step(_fe_problem.timeStep()), 46 57042 : _dt(_fe_problem.dt()), 47 57042 : _dt_old(_fe_problem.dtOld()), 48 57042 : _n_nonlinear_iterations(0), 49 57042 : _n_linear_iterations(0), 50 57042 : _is_lumped(false), 51 57042 : _var_restriction(declareRestartableData<bool>( 52 114084 : "var_restriction", !getParam<std::vector<VariableName>>("variables").empty())), 53 57042 : _local_indices(declareRestartableData<std::vector<dof_id_type>>("local_indices")), 54 57042 : _vars(declareRestartableData<std::unordered_set<unsigned int>>("vars")), 55 342252 : _from_subvector(NumericVector<Number>::build(this->comm())) 56 : { 57 57042 : _fe_problem.setUDotRequested(true); 58 57042 : } 59 : 60 : void 61 25482 : TimeIntegrator::init() 62 : { 63 25482 : if (!_var_restriction) 64 25370 : return; 65 : 66 112 : const auto & var_names = getParam<std::vector<VariableName>>("variables"); 67 112 : std::vector<unsigned int> var_num_vec; 68 112 : auto & lm_sys = _sys.system(); 69 112 : lm_sys.get_all_variable_numbers(var_num_vec); 70 112 : std::unordered_set<unsigned int> var_nums(var_num_vec.begin(), var_num_vec.end()); 71 224 : for (const auto & var_name : var_names) 72 112 : if (lm_sys.has_variable(var_name)) 73 : { 74 112 : const auto var_num = lm_sys.variable_number(var_name); 75 112 : _vars.insert(var_num); 76 112 : var_nums.erase(var_num); 77 : } 78 : 79 : // If var_nums is empty then that means the user has specified all the variables in this system 80 112 : if (var_nums.empty()) 81 : { 82 0 : _var_restriction = false; 83 0 : return; 84 : } 85 : 86 112 : std::vector<dof_id_type> var_dof_indices, work_vec; 87 224 : for (const auto var_num : _vars) 88 : { 89 112 : work_vec = _local_indices; 90 112 : _local_indices.clear(); 91 112 : lm_sys.get_dof_map().local_variable_indices(var_dof_indices, lm_sys.get_mesh(), var_num); 92 112 : std::merge(work_vec.begin(), 93 : work_vec.end(), 94 : var_dof_indices.begin(), 95 : var_dof_indices.end(), 96 : std::back_inserter(_local_indices)); 97 : } 98 : 99 112 : _solution_sub = NumericVector<Number>::build(_solution->comm()); 100 112 : _solution_old_sub = NumericVector<Number>::build(_solution_old.comm()); 101 112 : } 102 : 103 : void 104 0 : TimeIntegrator::solve() 105 : { 106 0 : mooseError("Calling TimeIntegrator::solve() is no longer supported"); 107 : } 108 : 109 : void 110 245224 : TimeIntegrator::setNumIterationsLastSolve() 111 : { 112 245224 : _n_nonlinear_iterations = getNumNonlinearIterationsLastSolve(); 113 245224 : _n_linear_iterations = getNumLinearIterationsLastSolve(); 114 245224 : } 115 : 116 : unsigned int 117 264095 : TimeIntegrator::getNumNonlinearIterationsLastSolve() const 118 : { 119 264095 : return _nonlinear_implicit_system->n_nonlinear_iterations(); 120 : } 121 : 122 : unsigned int 123 264095 : TimeIntegrator::getNumLinearIterationsLastSolve() const 124 : { 125 264095 : auto & nonlinear_solver = _nonlinear_implicit_system->nonlinear_solver; 126 : libmesh_assert(nonlinear_solver); 127 : 128 264095 : return nonlinear_solver->get_total_linear_iterations(); 129 : } 130 : 131 : void 132 9191 : TimeIntegrator::copyVector(const NumericVector<Number> & from, NumericVector<Number> & to) 133 : { 134 9191 : if (!_var_restriction) 135 8421 : to = from; 136 : else 137 : { 138 770 : auto to_sub = to.get_subvector(_local_indices); 139 770 : from.create_subvector(*_from_subvector, _local_indices, false); 140 770 : *to_sub = *_from_subvector; 141 770 : to.restore_subvector(std::move(to_sub), _local_indices); 142 770 : } 143 9191 : } 144 : 145 : bool 146 9356175 : TimeIntegrator::integratesVar(const unsigned int var_num) const 147 : { 148 9356175 : if (!_var_restriction) 149 9287991 : return true; 150 : 151 68184 : return _vars.count(var_num); 152 : } 153 : 154 : void 155 4423048 : TimeIntegrator::computeDuDotDu() 156 : { 157 4423048 : const auto coeff = duDotDuCoeff(); 158 13270226 : for (const auto i : index_range(_du_dot_du)) 159 8847178 : if (integratesVar(i)) 160 8830340 : _du_dot_du[i] = coeff / _dt; 161 4423048 : }