Line data Source code
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" 11 : #include "ElectromagneticConstants.h" 12 : #include <complex> 13 : 14 : registerMooseObject("ElectromagneticsApp", ReflectionCoefficient); 15 : 16 : InputParameters 17 46 : ReflectionCoefficient::validParams() 18 : { 19 46 : InputParameters params = SidePostprocessor::validParams(); 20 46 : 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 92 : params.addRequiredCoupledVar( 24 : "field_real", "The name of the real field variable this postprocessor operates on."); 25 92 : params.addRequiredCoupledVar("field_imag", "Coupled imaginary field variable."); 26 92 : params.addRequiredRangeCheckedParam<Real>("theta", "theta>=0", "Wave incidence angle"); 27 92 : params.addRequiredRangeCheckedParam<Real>("length", "length>0", "Domain length"); 28 92 : params.addRequiredRangeCheckedParam<Real>("k", "k>0", "Wave number"); 29 138 : params.addRangeCheckedParam<Real>( 30 92 : "incoming_field_magnitude", 1.0, "incoming_field_magnitude>0", "Incoming field magnitude"); 31 46 : return params; 32 0 : } 33 : 34 25 : ReflectionCoefficient::ReflectionCoefficient(const InputParameters & parameters) 35 : : SidePostprocessor(parameters), 36 : MooseVariableInterface<Real>(this, false, "field_real"), 37 25 : _qp(0), 38 25 : _coupled_real(coupledValue("field_real")), 39 25 : _coupled_imag(coupledValue("field_imag")), 40 50 : _theta(getParam<Real>("theta")), 41 50 : _length(getParam<Real>("length")), 42 50 : _k(getParam<Real>("k")), 43 100 : _incoming_mag(getParam<Real>("incoming_field_magnitude")) 44 : { 45 25 : } 46 : 47 : void 48 19 : ReflectionCoefficient::initialize() 49 : { 50 19 : if (_current_elem->dim() > 1) 51 : { 52 2 : 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 : 56 17 : _reflection_coefficient = 0; 57 17 : } 58 : 59 : void 60 9 : ReflectionCoefficient::execute() 61 : { 62 9 : _reflection_coefficient = computeReflection(); 63 9 : } 64 : 65 : PostprocessorValue 66 14 : ReflectionCoefficient::getValue() const 67 : { 68 14 : return _reflection_coefficient; 69 : } 70 : 71 : void 72 3 : ReflectionCoefficient::threadJoin(const UserObject & y) 73 : { 74 : const auto & pps = static_cast<const ReflectionCoefficient &>(y); 75 3 : Real temp_rc = _reflection_coefficient; 76 3 : _reflection_coefficient = std::max(temp_rc, pps._reflection_coefficient); 77 3 : } 78 : 79 : void 80 14 : ReflectionCoefficient::finalize() 81 : { 82 14 : gatherMax(_reflection_coefficient); 83 14 : } 84 : 85 : Real 86 9 : ReflectionCoefficient::computeReflection() 87 : { 88 9 : std::complex<double> field(_coupled_real[_qp], _coupled_imag[_qp]); 89 : 90 : std::complex<double> incoming_wave = 91 9 : _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 9 : return std::pow(std::abs(reflection_coefficient_complex), 2); 98 : }