https://mooseframework.inl.gov
ADConductionCurrent.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 "ADConductionCurrent.h"
11 #include "ElectromagneticEnums.h"
13 #include <complex>
14 
15 registerMooseObject("ElectromagneticsApp", ADConductionCurrent);
16 
19 {
21  params.addClassDescription(
22  "Calculates the current source term in the Helmholtz wave equation using "
23  "the conduction formulation of the current.");
24  params.addRequiredCoupledVar("E_real", "The real component of the electric field.");
25  params.deprecateParam("E_real", "field_real", "10/01/2025");
26  params.addRequiredCoupledVar("field_real", "The real component of the electric field.");
27  params.addRequiredCoupledVar("E_imag", "The imaginary component of the electric field.");
28  params.deprecateParam("E_imag", "field_imag", "10/01/2025");
29  params.addRequiredCoupledVar("field_imag", "The imaginary component of the electric field.");
30 
31  params.addParam<MaterialPropertyName>(
32  "conductivity_real", 1.0, "The real component of the material conductivity.");
33  params.addParam<MaterialPropertyName>(
34  "conductivity_imag", 0.0, "The imaginary component of the material conductivity.");
35 
36  params.addParam<MaterialPropertyName>(
37  "ang_freq_real", "ang_freq", "The real component of the angular drive frequency.");
38  params.addParam<MaterialPropertyName>(
39  "ang_freq_imag", 0.0, "The imaginary component of the angular drive frequency.");
40 
41  params.addParam<MaterialPropertyName>(
42  "permeability_real", "mu_vacuum", "The real component of the material permeability.");
43  params.addParam<MaterialPropertyName>(
44  "permeability_imag", 0.0, "The imaginary component of the material permeability.");
45 
46  MooseEnum component("real imaginary");
47  params.addParam<MooseEnum>("component", component, "Component of field (real or imaginary).");
48  return params;
49 }
50 
52  : ADVectorKernel(parameters),
53  _field_real(adCoupledVectorValue("field_real")),
54  _field_imag(adCoupledVectorValue("field_imag")),
55 
56  _cond_real(getADMaterialProperty<Real>("conductivity_real")),
57  _cond_imag(getADMaterialProperty<Real>("conductivity_imag")),
58 
59  _omega_real(getADMaterialProperty<Real>("ang_freq_real")),
60  _omega_imag(getADMaterialProperty<Real>("ang_freq_imag")),
61 
62  _mu_real(getADMaterialProperty<Real>("permeability_real")),
63  _mu_imag(getADMaterialProperty<Real>("permeability_imag")),
64 
65  _component(getParam<MooseEnum>("component"))
66 {
67 }
68 
69 ADReal
71 {
72  // TODO: In the future, need to add an AD capability to std::complex
73  ADReal mu_omega_real = _mu_real[_qp] * _omega_real[_qp] - _mu_imag[_qp] * _omega_imag[_qp];
74  ADReal mu_omega_imag = _mu_real[_qp] * _omega_imag[_qp] + _mu_imag[_qp] * _omega_real[_qp];
75 
76  ADRealVectorValue sigma_field_real =
78  ADRealVectorValue sigma_field_imag =
80 
81  if (_component == EM::REAL)
82  {
83  return _test[_i][_qp] * -1.0 *
84  (mu_omega_imag * sigma_field_real + mu_omega_real * sigma_field_imag);
85  }
86  else
87  {
88  return _test[_i][_qp] * (mu_omega_real * sigma_field_real - mu_omega_imag * sigma_field_imag);
89  }
90 }
Calculates the current source term in the Helmholtz wave equation using the conduction formulation of...
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
const ADMaterialProperty< Real > & _cond_imag
Imaginary component of the material conductivity (in S/m)
const ADVectorVariableValue & _field_imag
Vector variable of the imaginary component of the electric field.
static InputParameters validParams()
MooseEnum _component
Component of the field vector (real or imaginary)
const ADMaterialProperty< Real > & _omega_real
Real component of the angular drive frequency (in rad/s)
const ADTemplateVariableTestValue< T > & _test
DualNumber< Real, DNDerivativeType, true > ADReal
registerMooseObject("ElectromagneticsApp", ADConductionCurrent)
const ADVectorVariableValue & _field_real
Vector variable of the real component of the electric field.
void deprecateParam(const std::string &old_name, const std::string &new_name, const std::string &removal_date)
ADConductionCurrent(const InputParameters &parameters)
virtual ADReal computeQpResidual() override
unsigned int _i
const ADMaterialProperty< Real > & _cond_real
Real component of the material conductivity (in S/m)
const ADMaterialProperty< Real > & _omega_imag
Imaginary component of the angular drive frequency (in rad/s)
static InputParameters validParams()
void addRequiredCoupledVar(const std::string &name, const std::string &doc_string)
const ADMaterialProperty< Real > & _mu_real
Real component of the material permeability (in H/m)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void addClassDescription(const std::string &doc_string)
unsigned int _qp
const ADMaterialProperty< Real > & _mu_imag
Imaginary component of the material permeability (in H/m)