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 "MassLumpedTimeDerivative.h" 11 : 12 : // MOOSE includes 13 : #include "Assembly.h" 14 : #include "MooseVariableFE.h" 15 : 16 : #include "libmesh/quadrature.h" 17 : 18 : registerMooseObject("MooseApp", MassLumpedTimeDerivative); 19 : 20 : InputParameters 21 14542 : MassLumpedTimeDerivative::validParams() 22 : { 23 14542 : InputParameters params = TimeKernel::validParams(); 24 14542 : params.addClassDescription( 25 : "Lumped formulation of the time derivative $\\frac{\\partial u}{\\partial t}$. Its " 26 : "corresponding weak form is $\\dot{u_i}(\\psi_i, 1)$ where $\\dot{u_i}$ denotes the time " 27 : "derivative of the solution coefficient associated with node $i$."); 28 14542 : return params; 29 0 : } 30 : 31 144 : MassLumpedTimeDerivative::MassLumpedTimeDerivative(const InputParameters & parameters) 32 144 : : TimeKernel(parameters), _u_dot_nodal(_var.dofValuesDot()) 33 : { 34 144 : } 35 : 36 : Real 37 2133108 : MassLumpedTimeDerivative::computeQpResidual() 38 : { 39 2133108 : return _test[_i][_qp] * _u_dot_nodal[_i]; 40 : } 41 : 42 : Real 43 2132708 : MassLumpedTimeDerivative::computeQpJacobian() 44 : { 45 2132708 : return _test[_i][_qp] * _du_dot_du[_qp]; 46 : } 47 : 48 : void 49 140618 : MassLumpedTimeDerivative::computeJacobian() 50 : { 51 140618 : prepareMatrixTag(_assembly, _var.number(), _var.number()); 52 683506 : for (_i = 0; _i < _test.size(); _i++) 53 2675596 : for (_qp = 0; _qp < _qrule->n_points(); _qp++) 54 2132708 : _local_ke(_i, _i) += _JxW[_qp] * _coord[_qp] * computeQpJacobian(); 55 140618 : accumulateTaggedLocalMatrix(); 56 140618 : }