Skip to content

Commit

Permalink
Create hybrid CG-DG scheme
Browse files Browse the repository at this point in the history
  • Loading branch information
lindsayad committed Apr 13, 2023
1 parent 00f405f commit 3d896ad
Show file tree
Hide file tree
Showing 8 changed files with 287 additions and 0 deletions.
26 changes: 26 additions & 0 deletions framework/include/dgkernels/ADDGConvection.h
@@ -0,0 +1,26 @@
//* This file is part of the MOOSE framework
//* https://www.mooseframework.org
//*
//* All rights reserved, see COPYRIGHT for full restrictions
//* https://github.com/idaholab/moose/blob/master/COPYRIGHT
//*
//* Licensed under LGPL 2.1, please see LICENSE for details
//* https://www.gnu.org/licenses/lgpl-2.1.html

#pragma once

#include "ADDGKernel.h"

class ADDGConvection : public ADDGKernel
{
public:
static InputParameters validParams();

ADDGConvection(const InputParameters & parameters);

protected:
virtual ADReal computeQpResidual(Moose::DGResidualType type);

const ADMaterialProperty<RealVectorValue> & _velocity;
const ADMaterialProperty<RealVectorValue> & _velocity_neighbor;
};
30 changes: 30 additions & 0 deletions framework/include/kernels/ADConservativeAdvection.h
@@ -0,0 +1,30 @@
//* This file is part of the MOOSE framework
//* https://www.mooseframework.org
//*
//* All rights reserved, see COPYRIGHT for full restrictions
//* https://github.com/idaholab/moose/blob/master/COPYRIGHT
//*
//* Licensed under LGPL 2.1, please see LICENSE for details
//* https://www.gnu.org/licenses/lgpl-2.1.html

#pragma once

#include "ADKernel.h"

/**
* Advection of the variable by the velocity provided by the user.
* Options for numerical stabilization are: none; full upwinding
*/
class ADConservativeAdvection : public ADKernel
{
public:
static InputParameters validParams();

ADConservativeAdvection(const InputParameters & parameters);

protected:
virtual ADReal computeQpResidual() override;

/// advection velocity
const ADMaterialProperty<RealVectorValue> & _velocity;
};
52 changes: 52 additions & 0 deletions framework/src/dgkernels/ADDGConvection.C
@@ -0,0 +1,52 @@
//* This file is part of the MOOSE framework
//* https://www.mooseframework.org
//*
//* All rights reserved, see COPYRIGHT for full restrictions
//* https://github.com/idaholab/moose/blob/master/COPYRIGHT
//*
//* Licensed under LGPL 2.1, please see LICENSE for details
//* https://www.gnu.org/licenses/lgpl-2.1.html

#include "ADDGConvection.h"

registerMooseObject("MooseApp", ADDGConvection);

InputParameters
ADDGConvection::validParams()
{
InputParameters params = ADDGKernel::validParams();
params.addRequiredParam<MaterialPropertyName>("velocity", "Velocity vector");
params.addClassDescription("DG for convection");
return params;
}

ADDGConvection::ADDGConvection(const InputParameters & parameters)
: ADDGKernel(parameters),
_velocity(getADMaterialProperty<RealVectorValue>("velocity")),
_velocity_neighbor(getNeighborADMaterialProperty<RealVectorValue>("velocity"))
{
}

ADReal
ADDGConvection::computeQpResidual(Moose::DGResidualType type)
{
ADReal r = 0;

auto average = [](const auto & elem_value, const auto & neighbor_value)
{ return (elem_value + neighbor_value) / 2; };

const auto vdotn = average(_velocity[_qp], _velocity_neighbor[_qp]) * _normals[_qp];

switch (type)
{
case Moose::Element:
r += vdotn * average(_u[_qp], _u_neighbor[_qp]) * _test[_i][_qp];
break;

case Moose::Neighbor:
r += vdotn * average(_u[_qp], _u_neighbor[_qp]) * _test_neighbor[_i][_qp];
break;
}

return r;
}
34 changes: 34 additions & 0 deletions framework/src/kernels/ADConservativeAdvection.C
@@ -0,0 +1,34 @@
//* This file is part of the MOOSE framework
//* https://www.mooseframework.org
//*
//* All rights reserved, see COPYRIGHT for full restrictions
//* https://github.com/idaholab/moose/blob/master/COPYRIGHT
//*
//* Licensed under LGPL 2.1, please see LICENSE for details
//* https://www.gnu.org/licenses/lgpl-2.1.html

#include "ADConservativeAdvection.h"
#include "SystemBase.h"

registerMooseObject("MooseApp", ADConservativeAdvection);

InputParameters
ADConservativeAdvection::validParams()
{
InputParameters params = ADKernel::validParams();
params.addClassDescription("Conservative form of $\\nabla \\cdot \\vec{v} u$ which in its weak "
"form is given by: $(-\\nabla \\psi_i, \\vec{v} u)$.");
params.addRequiredParam<MaterialPropertyName>("velocity", "Velocity vector");
return params;
}

ADConservativeAdvection::ADConservativeAdvection(const InputParameters & parameters)
: ADKernel(parameters), _velocity(getADMaterialProperty<RealVectorValue>("velocity"))
{
}

ADReal
ADConservativeAdvection::computeQpResidual()
{
return -_grad_test[_i][_qp] * _velocity[_qp] * _u[_qp];
}
36 changes: 36 additions & 0 deletions modules/navier_stokes/include/kernels/CGMass.h
@@ -0,0 +1,36 @@
//* This file is part of the MOOSE framework
//* https://www.mooseframework.org
//*
//* All rights reserved, see COPYRIGHT for full restrictions
//* https://github.com/idaholab/moose/blob/master/COPYRIGHT
//*
//* Licensed under LGPL 2.1, please see LICENSE for details
//* https://www.gnu.org/licenses/lgpl-2.1.html

#pragma once

#include "ADKernel.h"

// Forward Declarations

/**
* This class computes the mass equation residual and Jacobian
* contributions for the incompressible Navier-Stokes momentum
* equation.
*/
class CGMass : public ADKernel
{
public:
static InputParameters validParams();

CGMass(const InputParameters & parameters);

virtual ~CGMass() {}

protected:
virtual ADReal computeQpResidual() override;

const ADVariableGradient & _grad_u_vel;
const ADVariableGradient & _grad_v_vel;
const ADVariableGradient & _grad_w_vel;
};
34 changes: 34 additions & 0 deletions modules/navier_stokes/include/kernels/DGMomentumPressure.h
@@ -0,0 +1,34 @@
//* This file is part of the MOOSE framework
//* https://www.mooseframework.org
//*
//* All rights reserved, see COPYRIGHT for full restrictions
//* https://github.com/idaholab/moose/blob/master/COPYRIGHT
//*
//* Licensed under LGPL 2.1, please see LICENSE for details
//* https://www.gnu.org/licenses/lgpl-2.1.html

#pragma once

#include "ADKernel.h"

// Forward Declarations

/**
* This class computes the mass equation residual and Jacobian
* contributions for the incompressible Navier-Stokes momentum
* equation.
*/
class DGMomentumPressure : public ADKernel
{
public:
static InputParameters validParams();

DGMomentumPressure(const InputParameters & parameters);

virtual ~DGMomentumPressure() {}

protected:
virtual ADReal computeQpResidual() override;

const ADVariableGradient & _grad_pressure;
};
40 changes: 40 additions & 0 deletions modules/navier_stokes/src/kernels/CGMass.C
@@ -0,0 +1,40 @@
//* This file is part of the MOOSE framework
//* https://www.mooseframework.org
//*
//* All rights reserved, see COPYRIGHT for full restrictions
//* https://github.com/idaholab/moose/blob/master/COPYRIGHT
//*
//* Licensed under LGPL 2.1, please see LICENSE for details
//* https://www.gnu.org/licenses/lgpl-2.1.html

#include "CGMass.h"
#include "Function.h"

registerMooseObject("NavierStokesApp", CGMass);

InputParameters
CGMass::validParams()
{
InputParameters params = ADKernel::validParams();
params.addRequiredCoupledVar("u", "x-velocity");
params.addCoupledVar("v", 0, "y-velocity"); // only required in 2D and 3D
params.addCoupledVar("w", 0, "z-velocity"); // only required in 3D
return params;
}

CGMass::CGMass(const InputParameters & parameters)
: ADKernel(parameters),
_grad_u_vel(adCoupledGradient("u")),
_grad_v_vel(isCoupled("v") ? adCoupledGradient("v") : _ad_grad_zero),
_grad_w_vel(isCoupled("w") ? adCoupledGradient("w") : _ad_grad_zero)
{
}

ADReal
CGMass::computeQpResidual()
{
// (div u) * q
// Note: we (arbitrarily) multiply this term by -1 so that it matches the -p(div v)
// term in the momentum equation. Not sure if that is really important?
return -(_grad_u_vel[_qp](0) + _grad_v_vel[_qp](1) + _grad_w_vel[_qp](2)) * _test[_i][_qp];
}
35 changes: 35 additions & 0 deletions modules/navier_stokes/src/kernels/DGMomentumPressure.C
@@ -0,0 +1,35 @@
//* This file is part of the MOOSE framework
//* https://www.mooseframework.org
//*
//* All rights reserved, see COPYRIGHT for full restrictions
//* https://github.com/idaholab/moose/blob/master/COPYRIGHT
//*
//* Licensed under LGPL 2.1, please see LICENSE for details
//* https://www.gnu.org/licenses/lgpl-2.1.html

#include "DGMomentumPressure.h"
#include "Function.h"

registerMooseObject("NavierStokesApp", DGMomentumPressure);

InputParameters
DGMomentumPressure::validParams()
{
InputParameters params = ADKernel::validParams();
params.addRequiredCoupledVar(NS::pressure, "The pressure variable");
params.addRequiredParam<unsigned short>("component", "The velocity component.");
return params;
}

DGMomentumPressure::DGMomentumPressure(const InputParameters & parameters)
: ADKernel(parameters),
_grad_pressure(adCoupledGradient(NS::pressure)),
_component(getParam<unsigned short>("component"))
{
}

ADReal
DGMomentumPressure::computeQpResidual()
{
return _grad_pressure[_qp](_component) * _test[_i][_qp];
}

0 comments on commit 3d896ad

Please sign in to comment.