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 "MatCoupledForce.h" 11 : 12 : #include "MooseVariable.h" 13 : 14 : registerMooseObject("MooseApp", MatCoupledForce); 15 : 16 : InputParameters 17 14290 : MatCoupledForce::validParams() 18 : { 19 14290 : InputParameters params = Kernel::validParams(); 20 : 21 14290 : params.addClassDescription( 22 : "Implements a forcing term RHS of the form PDE = RHS, where RHS = Sum_j c_j * m_j * v_j. " 23 : "c_j, m_j, and v_j are provided as real coefficients, material properties, and coupled " 24 : "variables, respectively."); 25 14290 : params.addRequiredCoupledVar("v", "The coupled variables which provide the force"); 26 14290 : params.addParam<std::vector<Real>>( 27 : "coef", "Coefficents ($\\sigma$) multiplier for the coupled force term."); 28 14290 : params.addParam<std::vector<MaterialPropertyName>>("material_properties", 29 : "The coupled material properties."); 30 14290 : return params; 31 0 : } 32 : 33 13 : MatCoupledForce::MatCoupledForce(const InputParameters & parameters) 34 : : Kernel(parameters), 35 13 : _n_coupled(coupledComponents("v")), 36 13 : _coupled_props(isParamValid("material_properties")), 37 13 : _v_var(coupledIndices("v")), 38 13 : _v(coupledValues("v")), 39 13 : _coef(isParamValid("coef") ? getParam<std::vector<Real>>("coef") 40 26 : : std::vector<Real>(_n_coupled, 1)) 41 : { 42 39 : for (MooseIndex(_n_coupled) j = 0; j < _n_coupled; ++j) 43 : { 44 26 : _v_var_to_index[_v_var[j]] = j; 45 : 46 26 : if (_var.number() == _v_var[j]) 47 0 : paramError("v", 48 : "Coupled variable 'v' needs to be different from 'variable' with MatCoupledForce, " 49 : "consider using Reaction or somethig similar"); 50 : } 51 : 52 13 : if (isParamValid("coef") && _coef.size() != _n_coupled) 53 0 : paramError("coef", "Size of coef must be equal to size of v"); 54 : 55 13 : if (_coupled_props) 56 : { 57 13 : _mat_props.resize(_n_coupled); 58 : std::vector<MaterialPropertyName> names = 59 13 : getParam<std::vector<MaterialPropertyName>>("material_properties"); 60 13 : if (names.size() != _n_coupled) 61 0 : paramError("material_properties", "Size must be equal to number of coupled variables"); 62 39 : for (unsigned int j = 0; j < _n_coupled; ++j) 63 26 : _mat_props[j] = &getMaterialPropertyByName<Real>(names[j]); 64 13 : } 65 13 : } 66 : 67 : Real 68 176000 : MatCoupledForce::computeQpResidual() 69 : { 70 176000 : Real r = 0; 71 176000 : if (_coupled_props) 72 528000 : for (unsigned int j = 0; j < _n_coupled; ++j) 73 352000 : r += -_coef[j] * (*_mat_props[j])[_qp] * (*_v[j])[_qp]; 74 : else 75 0 : for (unsigned int j = 0; j < _n_coupled; ++j) 76 0 : r += -_coef[j] * (*_v[j])[_qp]; 77 176000 : return r * _test[_i][_qp]; 78 : } 79 : 80 : Real 81 102400 : MatCoupledForce::computeQpJacobian() 82 : { 83 102400 : return 0; 84 : } 85 : 86 : Real 87 0 : MatCoupledForce::computeQpOffDiagJacobian(unsigned int jvar) 88 : { 89 0 : auto p = _v_var_to_index.find(jvar); 90 0 : if (p == _v_var_to_index.end()) 91 0 : return 0; 92 : 93 0 : if (_coupled_props) 94 0 : return -_coef[p->second] * (*_mat_props[p->second])[_qp] * _phi[_j][_qp] * _test[_i][_qp]; 95 0 : return -_coef[p->second] * _phi[_j][_qp] * _test[_i][_qp]; 96 : }