21 "CURRENTLY ONLY FOR 1D PLANE WAVE SOLVES. Calculate power reflection coefficient " 22 "for impinging wave on a surface. Assumes that wave of form F = F_incoming + R*F_reflected");
24 "field_real",
"The name of the real field variable this postprocessor operates on.");
30 "incoming_field_magnitude", 1.0,
"incoming_field_magnitude>0",
"Incoming field magnitude");
38 _coupled_real(coupledValue(
"field_real")),
39 _coupled_imag(coupledValue(
"field_imag")),
40 _theta(getParam<
Real>(
"theta")),
41 _length(getParam<
Real>(
"length")),
42 _k(getParam<
Real>(
"k")),
43 _incoming_mag(getParam<
Real>(
"incoming_field_magnitude"))
52 mooseError(
"The ReflectionCoefficient object is not currently configured to work with 2D or 3D " 53 "meshes. Please disable this object or setup a 1D mesh!");
90 std::complex<double> incoming_wave =
92 std::complex<double> reversed_wave =
95 std::complex<double> reflection_coefficient_complex = (field - incoming_wave) / reversed_wave;
97 return std::pow(std::abs(reflection_coefficient_complex), 2);
virtual Real computeReflection()
compute reflection coefficient
unsigned int _qp
quadrature point
const Real _theta
Wave incidence angle.
virtual PostprocessorValue getValue() const override
Real _reflection_coefficient
Value of the reflection coefficient.
const std::vector< double > y
const Real _k
Wave number.
const VariableValue & _coupled_imag
Imaginary component of the coupled field variable.
const Real _incoming_mag
Incoming field magnitude.
CURRENTLY ONLY FOR 1D PLANE WAVE SOLVES.
static InputParameters validParams()
virtual void execute() override
virtual void initialize() override
const Real _length
Domain length.
virtual void threadJoin(const UserObject &y) override
virtual void finalize() override
registerMooseObject("ElectromagneticsApp", ReflectionCoefficient)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static InputParameters validParams()
void mooseError(Args &&... args) const
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
ReflectionCoefficient(const InputParameters ¶meters)
const Elem *const & _current_elem
MooseUnits pow(const MooseUnits &, int)
const VariableValue & _coupled_real
Real component of the coupled field variable.