https://mooseframework.inl.gov
VectorCurrentSource.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 "VectorCurrentSource.h"
11 #include "ElectromagneticEnums.h"
13 #include "Function.h"
14 #include <complex>
15 
16 registerMooseObject("ElectromagneticsApp", VectorCurrentSource);
17 
20 {
22  params.addClassDescription(
23  "Kernel to calculate the current source term in the Helmholtz wave equation.");
24  params.addParam<FunctionName>("function_coefficient",
25  1.0,
26  "Function coefficient multiplier for current source (normally "
27  "$\\omega$ or $\\omega \\cdot \\mu$).");
28  params.addRequiredParam<FunctionName>("source_real", "Current Source vector, real component");
29  params.addRequiredParam<FunctionName>("source_imag",
30  "Current Source vector, imaginary component");
31  MooseEnum component("real imaginary");
32  params.addParam<MooseEnum>("component", component, "Component of field (real or imaginary).");
33  return params;
34 }
35 
37  : VectorKernel(parameters),
38  _func(getFunction("function_coefficient")),
39  _source_real(getFunction("source_real")),
40  _source_imag(getFunction("source_imag")),
41  _component(getParam<MooseEnum>("component"))
42 {
43 }
44 
45 Real
47 {
48  std::complex<double> source_0(_source_real.vectorValue(_t, _q_point[_qp])(0),
50  std::complex<double> source_1(_source_real.vectorValue(_t, _q_point[_qp])(1),
52  std::complex<double> source_2(_source_real.vectorValue(_t, _q_point[_qp])(2),
54  VectorValue<std::complex<double>> source(source_0, source_1, source_2);
55 
56  std::complex<double> res = EM::j * _func.value(_t, _q_point[_qp]) * source * _test[_i][_qp];
57 
58  if (_component == EM::REAL)
59  return res.real();
60  else
61  return res.imag();
62 }
63 
64 Real
66 {
67  return 0.0;
68 }
virtual Real computeQpJacobian() override
const Function & _func
Function coefficient.
const VectorVariableTestValue & _test
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
static const std::string component
Definition: NS.h:153
static InputParameters validParams()
const Function & _source_real
Real component of the current source.
void addRequiredParam(const std::string &name, const std::string &doc_string)
virtual Real computeQpResidual() override
VectorCurrentSource(const InputParameters &parameters)
unsigned int _i
MooseEnum _component
Component of the field vector (real or imaginary)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Function & _source_imag
Imaginary component of the current source.
static InputParameters validParams()
registerMooseObject("ElectromagneticsApp", VectorCurrentSource)
virtual RealVectorValue vectorValue(Real t, const Point &p) const
void addClassDescription(const std::string &doc_string)
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
Calculates the current source term contribution in the Helmholtz wave equation.
virtual Real value(Real t, const Point &p) const
const MooseArray< Point > & _q_point
unsigned int _qp