https://mooseframework.inl.gov
ReflectionCoefficient.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 "ReflectionCoefficient.h"
12 #include <complex>
13 
14 registerMooseObject("ElectromagneticsApp", ReflectionCoefficient);
15 
18 {
20  params.addClassDescription(
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");
23  params.addRequiredCoupledVar(
24  "field_real", "The name of the real field variable this postprocessor operates on.");
25  params.addRequiredCoupledVar("field_imag", "Coupled imaginary field variable.");
26  params.addRequiredRangeCheckedParam<Real>("theta", "theta>=0", "Wave incidence angle");
27  params.addRequiredRangeCheckedParam<Real>("length", "length>0", "Domain length");
28  params.addRequiredRangeCheckedParam<Real>("k", "k>0", "Wave number");
29  params.addRangeCheckedParam<Real>(
30  "incoming_field_magnitude", 1.0, "incoming_field_magnitude>0", "Incoming field magnitude");
31  return params;
32 }
33 
35  : SidePostprocessor(parameters),
36  MooseVariableInterface<Real>(this, false, "field_real"),
37  _qp(0),
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"))
44 {
45 }
46 
47 void
49 {
50  if (_current_elem->dim() > 1)
51  {
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!");
54  }
55 
57 }
58 
59 void
61 {
63 }
64 
67 {
69 }
70 
71 void
73 {
74  const auto & pps = static_cast<const ReflectionCoefficient &>(y);
75  Real temp_rc = _reflection_coefficient;
76  _reflection_coefficient = std::max(temp_rc, pps._reflection_coefficient);
77 }
78 
79 void
81 {
83 }
84 
85 Real
87 {
88  std::complex<double> field(_coupled_real[_qp], _coupled_imag[_qp]);
89 
90  std::complex<double> incoming_wave =
91  _incoming_mag * std::exp(EM::j * _k * _length * std::cos(_theta * libMesh::pi / 180.));
92  std::complex<double> reversed_wave =
93  _incoming_mag * std::exp(-EM::j * _k * _length * std::cos(_theta * libMesh::pi / 180.));
94 
95  std::complex<double> reflection_coefficient_complex = (field - incoming_wave) / reversed_wave;
96 
97  return std::pow(std::abs(reflection_coefficient_complex), 2);
98 }
virtual Real computeReflection()
compute reflection coefficient
unsigned int _qp
quadrature point
const Real _theta
Wave incidence angle.
void addRequiredRangeCheckedParam(const std::string &name, const std::string &parsed_function, const std::string &doc_string)
virtual PostprocessorValue getValue() const override
Real _reflection_coefficient
Value of the reflection coefficient.
void gatherMax(T &value)
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.
Real PostprocessorValue
virtual void threadJoin(const UserObject &y) override
virtual void finalize() override
registerMooseObject("ElectromagneticsApp", ReflectionCoefficient)
void addRequiredCoupledVar(const std::string &name, const std::string &doc_string)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static InputParameters validParams()
void mooseError(Args &&... args) const
void addClassDescription(const std::string &doc_string)
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
void addRangeCheckedParam(const std::string &name, const T &value, const std::string &parsed_function, const std::string &doc_string)
ReflectionCoefficient(const InputParameters &parameters)
const Elem *const & _current_elem
MooseUnits pow(const MooseUnits &, int)
const VariableValue & _coupled_real
Real component of the coupled field variable.
const Real pi