https://mooseframework.inl.gov
MFEMCurlAux.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 #ifdef MOOSE_MFEM_ENABLED
11 
12 #include "MFEMCurlAux.h"
13 #include "MFEMProblem.h"
14 
15 registerMooseObject("MooseApp", MFEMCurlAux);
16 
19 {
21  params.addClassDescription(
22  "Calculates the curl of an H(curl) conforming ND source variable and stores the result"
23  " on an H(div) conforming RT result auxvariable");
24  params.addRequiredParam<VariableName>("source",
25  "Vector H(curl) MFEMVariable to take the curl of.");
26  params.addParam<mfem::real_t>("scale_factor", 1.0, "Factor to scale result auxvariable by.");
27  return params;
28 }
29 
31  : MFEMAuxKernel(parameters),
32  _source_var_name(getParam<VariableName>("source")),
33  _source_var(*getMFEMProblem().getProblemData().gridfunctions.Get(_source_var_name)),
34  _scale_factor(getParam<mfem::real_t>("scale_factor")),
35  _curl(_source_var.ParFESpace(), _result_var.ParFESpace())
36 {
37  _curl.Assemble();
38  _curl.Finalize();
39 }
40 
41 // Computes the auxvariable.
42 void
44 {
45  _result_var = 0.0;
47 }
48 
49 #endif
const mfem::real_t _scale_factor
Scalar factor to multiply the result by.
Definition: MFEMCurlAux.h:39
MFEMCurlAux(const InputParameters &parameters)
Definition: MFEMCurlAux.C:30
Class to set an H(div) auxvariable to be the curl of an H(curl) vector variable.
Definition: MFEMCurlAux.h:21
virtual void execute() override
Computes the auxvariable.
Definition: MFEMCurlAux.C:43
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
registerMooseObject("MooseApp", MFEMCurlAux)
void addRequiredParam(const std::string &name, const std::string &doc_string)
This method adds a parameter and documentation string to the InputParameters object that will be extr...
const mfem::ParGridFunction & _source_var
Reference to source gridfunction.
Definition: MFEMCurlAux.h:37
mfem::ParGridFunction & _result_var
Reference to result gridfunction.
Definition: MFEMAuxKernel.h:36
Class to construct an auxiliary solver used to update an auxvariable.
Definition: MFEMAuxKernel.h:20
static InputParameters validParams()
Definition: MFEMAuxKernel.C:16
static InputParameters validParams()
Definition: MFEMCurlAux.C:18
mfem::common::ParDiscreteCurlOperator _curl
Curl operator.
Definition: MFEMCurlAux.h:41
void addClassDescription(const std::string &doc_string)
This method adds a description of the class that will be displayed in the input file syntax dump...
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an optional parameter and a documentation string to the InputParameters object...