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 : #ifdef MFEM_ENABLED 11 : 12 : #include "MFEMLinearElasticityKernel.h" 13 : #include "MFEMProblem.h" 14 : 15 : registerMooseObject("MooseApp", MFEMLinearElasticityKernel); 16 : 17 : InputParameters 18 8648 : MFEMLinearElasticityKernel::validParams() 19 : { 20 8648 : InputParameters params = MFEMKernel::validParams(); 21 8648 : params.addClassDescription( 22 : "The isotropic linear elasticity operator with weak form " 23 : "$(c_{ikjl} \\nabla u_j, \\nabla v_i)$, to be added to an MFEM problem, where " 24 : "$c_{ikjl}$ is the isotropic elasticity tensor, " 25 : "$c_{ikjl} = \\lambda \\delta_{ik} \\delta_{jl} + \\mu \\left( \\delta_{ij} \\delta_{kl} + " 26 : "\\delta_{il} \\delta_{jk} \\right)$, " 27 : "$\\lambda$ is the first Lame parameter, $\\lambda = \\frac{E\\nu}{(1-2\\nu)(1+\\nu)}$, " 28 : "$\\mu$ is the second Lame parameter, $\\mu = \\frac{E}{2(1+\\nu)}$, " 29 : "where $E$ is Young's modulus and $\\nu$ is Poisson's ratio."); 30 : 31 8648 : params.addParam<MFEMScalarCoefficientName>( 32 : "lambda", "1.", "Name of MFEM Lame constant lambda to multiply the div(u)*I term by"); 33 8648 : params.addParam<MFEMScalarCoefficientName>( 34 : "mu", "1.", "Name of MFEM Lame constant mu to multiply the gradients term by"); 35 : 36 8648 : return params; 37 0 : } 38 : 39 9 : MFEMLinearElasticityKernel::MFEMLinearElasticityKernel(const InputParameters & parameters) 40 : : MFEMKernel(parameters), 41 9 : _lambda_name(getParam<MFEMScalarCoefficientName>("lambda")), 42 9 : _mu_name(getParam<MFEMScalarCoefficientName>("mu")), 43 9 : _lambda(getScalarCoefficient(_lambda_name)), 44 18 : _mu(getScalarCoefficient(_mu_name)) 45 : { 46 9 : } 47 : 48 : mfem::BilinearFormIntegrator * 49 9 : MFEMLinearElasticityKernel::createBFIntegrator() 50 : { 51 9 : return new mfem::ElasticityIntegrator(_lambda, _mu); 52 : } 53 : 54 : #endif