https://mooseframework.inl.gov
JinSlabCoeffFunc.C
Go to the documentation of this file.
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 "JinSlabCoeffFunc.h"
11 #include "ElectromagneticEnums.h"
13 #include <complex>
14 
15 registerMooseObject("ElectromagneticsTestApp", JinSlabCoeffFunc);
16 
19 {
21  params.addClassDescription(
22  "Function describing a wave incident on a surface at a given angle, wavenumber, and domain "
23  "length, for use in the slab reflection benchmark.");
24  params.addRequiredParam<Real>("k", "Wavenumber");
25  params.addRequiredParam<Real>("theta", "Wave Incidence angle, in degrees");
26  params.addRequiredParam<Real>("length", "Length of slab domain");
27  params.addParam<Real>("coef", 1.0, "Real-valued function coefficient.");
28  MooseEnum component("real imaginary");
29  params.addParam<MooseEnum>("component", component, "Real or Imaginary wave component");
30  return params;
31 }
32 
34  : Function(parameters),
35  _k(getParam<Real>("k")),
36  _theta(getParam<Real>("theta")),
37  _length(getParam<Real>("length")),
38  _coef(getParam<Real>("coef")),
39  _component(getParam<MooseEnum>("component"))
40 {
41 }
42 
43 Real
44 JinSlabCoeffFunc::value(Real /*t*/, const Point & p) const
45 {
46  std::complex<double> eps_r = 4.0 + (2.0 - EM::j * 0.1) * std::pow((1 - p(0) / _length), 2);
47  std::complex<double> mu_r(2, -0.1);
48 
49  std::complex<double> val = _coef * mu_r * std::pow(_k, 2) *
50  (eps_r - (1.0 / mu_r) * std::pow(sin(_theta * libMesh::pi / 180.), 2));
51 
52  if (_component == EM::REAL)
53  return val.real();
54  else
55  return val.imag();
56 }
CTSub CT_OPERATOR_BINARY CTMul CTCompareLess CTCompareGreater CTCompareEqual _arg template * sin(_arg) *_arg.template D< dtag >()) CT_SIMPLE_UNARY_FUNCTION(tan
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
const Real _theta
Wave incidence angle, in degrees.
static const std::string component
Definition: NS.h:153
const Real _length
Length of slab domain (m)
const MooseEnum _component
Enum signifying the component of the function being calculated.
void addRequiredParam(const std::string &name, const std::string &doc_string)
Function for field coefficient in slab reflection benchmark case.
static InputParameters validParams()
JinSlabCoeffFunc(const InputParameters &parameters)
registerMooseObject("ElectromagneticsTestApp", JinSlabCoeffFunc)
virtual Real value(Real t, const Point &p) const override
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Real _coef
Real-valued function coefficient.
void addClassDescription(const std::string &doc_string)
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
static InputParameters validParams()
const Real _k
Wavenumber of incoming wave (1/m)
MooseUnits pow(const MooseUnits &, int)
const Real pi