MFEMComplexKernel

Summary

Base class for complex MFEM kernels applied to the weak form being solved.

Overview

Complex MFEM kernels are responsible for providing complex-valued domain integrators to add to the weak form of the FE problem accumulated in EquationSystem, along with any required marker arrays to restrict the integrator(s) to subdomains. The MFEMComplexKernel is a container for two objects of MFEMKernel type, one representing a real contribution to the equation system, and the other representing an imaginary contribution. The real and imaginary kernels need not be the same, but they do need to be applied to the same variable, which must be of MFEMComplexVariable type.

The integrators are applied to the weak form equation that is labeled according to the test variable name of the kernel returned from getTestVariableName(). In the case of bilinear (or nonlinear) forms, the trial variable that the integrator acts on is the variable returned from getTrialVariableName().

Example Input File Syntax

The real and imaginary contributions to the MFEMComplexKernel can be set up by using the sub-blocks RealComponent and ImagComponent on the script.

[Kernels<<<{"href": "../../../syntax/Kernels/index.html"}>>>]
  [diffusion_complex]
    type = MFEMComplexKernel<<<{"description": "Holds MFEMKernel objects for the real and imaginary parts of a complex kernel.", "href": "MFEMComplexKernel.html"}>>>
    variable<<<{"description": "Variable labelling the weak form this kernel is added to"}>>> = u
    [RealComponent]
      type = MFEMDiffusionKernel
      coefficient = stiffnessCoef
    []
    [ImagComponent]
      type = MFEMDiffusionKernel
      coefficient = 0.0
    []
  []
  [mass_complex]
    type = MFEMComplexKernel<<<{"description": "Holds MFEMKernel objects for the real and imaginary parts of a complex kernel.", "href": "MFEMComplexKernel.html"}>>>
    variable<<<{"description": "Variable labelling the weak form this kernel is added to"}>>> = u
    [RealComponent]
      type = MFEMMassKernel
      coefficient = massCoef
    []
    [ImagComponent]
      type = MFEMMassKernel
      coefficient = lossCoef
    []
  []
[]
(test/tests/mfem/complex/complex.i)

Note that the variable to which the kernel is applied must be complex-valued:

[Variables<<<{"href": "../../../syntax/Variables/index.html"}>>>]
  [u]
    type = MFEMComplexVariable<<<{"description": "Class for adding complex MFEM variables to the problem (`mfem::ParComplexGridFunction`s).", "href": "../variables/MFEMComplexVariable.html"}>>>
    fespace<<<{"description": "The finite element space this variable is defined on."}>>> = H1FESpace
  []
[]
(test/tests/mfem/complex/complex.i)

Input Parameters

  • blockThe list of subdomains (names or ids) that this object will be restricted to. Leave empty to apply to all subdomains.

    C++ Type:std::vector<SubdomainName>

    Controllable:No

    Description:The list of subdomains (names or ids) that this object will be restricted to. Leave empty to apply to all subdomains.

  • variableVariable labelling the weak form this kernel is added to

    C++ Type:VariableName

    Unit:(no unit assumed)

    Controllable:No

    Description:Variable labelling the weak form this kernel is added to

Optional Parameters

  • control_tagsAdds user-defined labels for accessing object parameters via control logic.

    C++ Type:std::vector<std::string>

    Controllable:No

    Description:Adds user-defined labels for accessing object parameters via control logic.

  • enableTrueSet the enabled status of the MooseObject.

    Default:True

    C++ Type:bool

    Controllable:No

    Description:Set the enabled status of the MooseObject.

Advanced Parameters

Input Files