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 "ExplicitMassDamping.h" 11 : #include "MooseError.h" 12 : #include "FEProblemBase.h" 13 : #include "NonlinearSystemBase.h" 14 : 15 : registerMooseObject("SolidMechanicsApp", ExplicitMassDamping); 16 : 17 : InputParameters 18 18 : ExplicitMassDamping::validParams() 19 : { 20 18 : InputParameters params = NodalKernel::validParams(); 21 18 : params.addClassDescription( 22 : "Adds Rayleigh mass damping, eta * M * v, to an ExplicitMixedOrder solid mechanics model."); 23 54 : params.addRangeCheckedParam<Real>( 24 : "eta", 25 36 : 0.0, 26 : "eta >= 0.0", 27 : "Damping strength. Mass damping helps to damp low-frequency, large-scale oscillations. A " 28 : "reasonable value to use initially is eta = 2 * d * omega, where d is the damping fraction " 29 : "(typically 0.1) and omega is the lowest relevant frequency mode of the model. For " 30 : "instance, if the model is can oscillate longitudinally, omega = (pi / model_length) * " 31 : "sqrt(youngs_modulus / density)"); 32 18 : return params; 33 0 : } 34 : 35 9 : ExplicitMassDamping::ExplicitMassDamping(const InputParameters & parameters) 36 : : NodalKernel(parameters), 37 9 : _eta(getParam<Real>("eta")), 38 9 : _mass_matrix_lumped(initLumpedMass()), 39 18 : _u_dot_old(_var.nodalValueDotOld()) 40 : { 41 9 : } 42 : 43 : Real 44 96 : ExplicitMassDamping::computeQpResidual() 45 : { 46 96 : const auto dofnum = _variable->nodalDofIndex(); // dof for current var 47 96 : return _eta * _mass_matrix_lumped(dofnum) * _u_dot_old; 48 : } 49 : 50 : const NumericVector<Number> & 51 9 : ExplicitMassDamping::initLumpedMass() 52 : { 53 9 : const auto & nl = _fe_problem.getNonlinearSystemBase(_sys.number()); 54 18 : if (nl.hasVector("mass_matrix_lumped")) 55 18 : return nl.getVector("mass_matrix_lumped"); 56 : 57 0 : mooseError("Lumped mass matrix is missing. Make sure ExplicitMixedOrder is being used as the " 58 : "time integrator."); 59 : }