Skip to content

Commit

Permalink
Add divergence and gradient auxkernels
Browse files Browse the repository at this point in the history
  • Loading branch information
alexanderianblair committed Dec 9, 2024
1 parent 2732f7b commit c0ccb4f
Show file tree
Hide file tree
Showing 6 changed files with 133 additions and 3 deletions.
4 changes: 2 additions & 2 deletions include/auxkernels/MFEMCurlAux.h
Original file line number Diff line number Diff line change
Expand Up @@ -3,7 +3,7 @@
#include "MFEMAuxKernel.h"

/*
Class to set an auxvariable to be the curl of a variable.
Class to set an H(div) auxvariable to be the curl of a H(curl) vector variable.
*/
class MFEMCurlAux : public MFEMAuxKernel
{
Expand All @@ -18,7 +18,7 @@ class MFEMCurlAux : public MFEMAuxKernel
virtual void execute() override;

protected:
// Name of auxvariable to store the result of the auxkernel in.
// Name of source MFEMVariable to take the curl of.
VariableName _source_var_name;
// Pointer to source gridfunction.
mfem::ParGridFunction & _source_var;
Expand Down
30 changes: 30 additions & 0 deletions include/auxkernels/MFEMDivAux.h
Original file line number Diff line number Diff line change
@@ -0,0 +1,30 @@
#pragma once
#include "../common/pfem_extras.hpp"
#include "MFEMAuxKernel.h"

/*
Class to set an L2 auxvariable to be the divergence of a H(div) vector variable.
*/
class MFEMDivAux : public MFEMAuxKernel
{
public:
static InputParameters validParams();

MFEMDivAux(const InputParameters & parameters);

virtual ~MFEMDivAux() = default;

// Computes the auxvariable.
virtual void execute() override;

protected:
// Name of source MFEMVariable to take the divergence of.
VariableName _source_var_name;
// Pointer to source gridfunction.
mfem::ParGridFunction & _source_var;
// FESpaces
mfem::ParFiniteElementSpace & _hdiv_fespace;
mfem::ParFiniteElementSpace & _l2_fespace;
// Divergence operator
mfem::common::ParDiscreteCurlOperator _div;
};
30 changes: 30 additions & 0 deletions include/auxkernels/MFEMGradAux.h
Original file line number Diff line number Diff line change
@@ -0,0 +1,30 @@
#pragma once
#include "../common/pfem_extras.hpp"
#include "MFEMAuxKernel.h"

/*
Class to set an H(curl) auxvariable to be the gradient of a H1 scalar variable.
*/
class MFEMGradAux : public MFEMAuxKernel
{
public:
static InputParameters validParams();

MFEMGradAux(const InputParameters & parameters);

virtual ~MFEMGradAux() = default;

// Computes the auxvariable.
virtual void execute() override;

protected:
// Name of source MFEMVariable to take the gradient of.
VariableName _source_var_name;
// Pointer to source gridfunction.
mfem::ParGridFunction & _source_var;
// FESpaces
mfem::ParFiniteElementSpace & _h1_fespace;
mfem::ParFiniteElementSpace & _hcurl_fespace;
// Grad operator
mfem::common::ParDiscreteGradOperator _grad;
};
2 changes: 1 addition & 1 deletion src/auxkernels/MFEMCurlAux.C
Original file line number Diff line number Diff line change
Expand Up @@ -4,7 +4,7 @@
registerMooseObject("PlatypusApp", MFEMCurlAux);

/*
Class to set an auxvariable to be the curl of a variable.
Class to set an H(div) auxvariable to be the curl of a H(curl) vector variable.
*/
InputParameters
MFEMCurlAux::validParams()
Expand Down
35 changes: 35 additions & 0 deletions src/auxkernels/MFEMDivAux.C
Original file line number Diff line number Diff line change
@@ -0,0 +1,35 @@
#pragma once
#include "MFEMDivAux.h"

registerMooseObject("PlatypusApp", MFEMDivAux);

/*
Class to set an L2 auxvariable to be the divergence of a H(div) vector variable.
*/
InputParameters
MFEMDivAux::validParams()
{
InputParameters params = MFEMAuxKernel::validParams();
params.addRequiredParam<VariableName>("source",
"Vector H(div) MFEMVariable to take the divergence of.");
return params;
}

MFEMDivAux::MFEMDivAux(const InputParameters & parameters)

Check warning on line 18 in src/auxkernels/MFEMDivAux.C

View check run for this annotation

Codecov / codecov/patch

src/auxkernels/MFEMDivAux.C#L18

Added line #L18 was not covered by tests
: MFEMAuxKernel(parameters),
_source_var_name(getParam<VariableName>("source")),
_source_var(*getMFEMProblem().getProblemData()._gridfunctions.Get(_source_var_name)),
_hdiv_fespace(*_source_var.ParFESpace()),
_l2_fespace(*_result_var.ParFESpace()),
_div(&_hdiv_fespace, &_l2_fespace)

Check warning on line 24 in src/auxkernels/MFEMDivAux.C

View check run for this annotation

Codecov / codecov/patch

src/auxkernels/MFEMDivAux.C#L20-L24

Added lines #L20 - L24 were not covered by tests
{
_div.Assemble();
_div.Finalize();

Check warning on line 27 in src/auxkernels/MFEMDivAux.C

View check run for this annotation

Codecov / codecov/patch

src/auxkernels/MFEMDivAux.C#L26-L27

Added lines #L26 - L27 were not covered by tests
}

// Computes the auxvariable.
void
MFEMDivAux::execute()

Check warning on line 32 in src/auxkernels/MFEMDivAux.C

View check run for this annotation

Codecov / codecov/patch

src/auxkernels/MFEMDivAux.C#L32

Added line #L32 was not covered by tests
{
_div.Mult(_source_var, _result_var);

Check warning on line 34 in src/auxkernels/MFEMDivAux.C

View check run for this annotation

Codecov / codecov/patch

src/auxkernels/MFEMDivAux.C#L34

Added line #L34 was not covered by tests
}
35 changes: 35 additions & 0 deletions src/auxkernels/MFEMGradAux.C
Original file line number Diff line number Diff line change
@@ -0,0 +1,35 @@
#pragma once
#include "MFEMGradAux.h"

registerMooseObject("PlatypusApp", MFEMGradAux);

/*
Class to set an H(curl) auxvariable to be the gradient of a H1 scalar variable.
*/
InputParameters
MFEMGradAux::validParams()
{
InputParameters params = MFEMAuxKernel::validParams();
params.addRequiredParam<VariableName>("source",
"Scalar H1 MFEMVariable to take the gradient of.");
return params;
}

MFEMGradAux::MFEMGradAux(const InputParameters & parameters)

Check warning on line 18 in src/auxkernels/MFEMGradAux.C

View check run for this annotation

Codecov / codecov/patch

src/auxkernels/MFEMGradAux.C#L18

Added line #L18 was not covered by tests
: MFEMAuxKernel(parameters),
_source_var_name(getParam<VariableName>("source")),
_source_var(*getMFEMProblem().getProblemData()._gridfunctions.Get(_source_var_name)),
_h1_fespace(*_source_var.ParFESpace()),
_hcurl_fespace(*_result_var.ParFESpace()),
_grad(&_h1_fespace, &_hcurl_fespace)

Check warning on line 24 in src/auxkernels/MFEMGradAux.C

View check run for this annotation

Codecov / codecov/patch

src/auxkernels/MFEMGradAux.C#L20-L24

Added lines #L20 - L24 were not covered by tests
{
_grad.Assemble();
_grad.Finalize();

Check warning on line 27 in src/auxkernels/MFEMGradAux.C

View check run for this annotation

Codecov / codecov/patch

src/auxkernels/MFEMGradAux.C#L26-L27

Added lines #L26 - L27 were not covered by tests
}

// Computes the auxvariable.
void
MFEMGradAux::execute()

Check warning on line 32 in src/auxkernels/MFEMGradAux.C

View check run for this annotation

Codecov / codecov/patch

src/auxkernels/MFEMGradAux.C#L32

Added line #L32 was not covered by tests
{
_grad.Mult(_source_var, _result_var);

Check warning on line 34 in src/auxkernels/MFEMGradAux.C

View check run for this annotation

Codecov / codecov/patch

src/auxkernels/MFEMGradAux.C#L34

Added line #L34 was not covered by tests
}

0 comments on commit c0ccb4f

Please sign in to comment.