add: enum for form of convection operator

This commit is contained in:
timovanopstal 2016-12-09 15:50:32 +01:00 committed by Arne Morten Kvarving
parent fccea559f1
commit 48801295b5
4 changed files with 32 additions and 15 deletions

View File

@ -34,7 +34,7 @@ void CompatibleOperators::Weak::Convection(std::vector<Matrix>& EM,
const Vec3& U,
const Tensor& dUdX,
double scale,
bool conservative)
WeakOperators::ConvectionForm form)
{
// Convection
static const double vidx[3][3] = {{1, 6, 7},
@ -121,7 +121,7 @@ void CompatibleOperators::Weak::Source(Vectors& EV,
void CompatibleOperators::Residual::Convection(Vectors& EV, const FiniteElement& fe,
const Vec3& U, const Tensor& dUdX,
const Vec3& UC, double scale,
bool conservative)
WeakOperators::ConvectionForm form)
{
size_t nsd = fe.grad(1).cols();
for (size_t n=1; n <= nsd; ++n) {

View File

@ -18,6 +18,7 @@ class Tensor;
class Vec3;
#include "MatVec.h"
#include "EqualOrderOperators.h"
/*! \brief Common operators using div-compatible discretizations.
* \details The operators use the block ordering used in the BlockElmMats class.
@ -41,11 +42,11 @@ public:
//! \param[out] EM The element matrix to add contribution to
//! \param[in] fe The finite element to evaluate for
//! \param[in] U Advecting field
//! \param[in] conservative True to use the conservative formulation
//! \param[in] form Which form of the convective form to use
//! \param[in] basis Basis to use
static void Convection(std::vector<Matrix>& EM, const FiniteElement& fe,
const Vec3& U, const Tensor& dUdX, double scale,
bool conservative=false);
WeakOperators::ConvectionForm form=WeakOperators::CONVECTIVE);
//! \brief Compute a gradient term.
//! \param[out] EM The element matrix to add contribution to
@ -104,11 +105,11 @@ public:
//! \param[in] dUdX Advected field gradient
//! \param[in] UC Advecting field
//! \param[in] scale Scaling factor for contribution
//! \param[in] conservative True to use the conservative formulation
//! \param[in] form Which form of the convective form to use
//! \param[in] basis Basis to use
static void Convection(Vectors& EV, const FiniteElement& fe,
const Vec3& U, const Tensor& dUdX, const Vec3& UC,
double scale, bool conservative=false);
const Vec3& U, const Tensor& dUdX, const Vec3& UC, double scale,
WeakOperators::ConvectionForm form=WeakOperators::CONVECTIVE);
//! \brief Compute a laplacian term in a residual vector.
//! \param[out] EV The element vector to add contribution to

View File

@ -72,10 +72,10 @@ void EqualOrderOperators::Weak::Advection(Matrix& EM, const FiniteElement& fe,
void EqualOrderOperators::Weak::Convection(Matrix& EM, const FiniteElement& fe,
const Vec3& U, const Tensor& dUdX,
double scale, bool conservative, int basis)
double scale, WeakOperators::ConvectionForm form, int basis)
{
size_t cmp = EM.rows() / fe.basis(basis).size();
if (conservative) {
if (form == WeakOperators::CONSERVATIVE) {
Advection(EM, fe, U, -scale, basis);
for (size_t i = 1;i <= fe.basis(basis).size();i++)
for (size_t j = 1;j <= fe.basis(basis).size();j++) {
@ -200,10 +200,10 @@ void EqualOrderOperators::Weak::Source(Vector& EV, const FiniteElement& fe,
void EqualOrderOperators::Residual::Convection(Vector& EV, const FiniteElement& fe,
const Vec3& U, const Tensor& dUdX,
const Vec3& UC, double scale,
bool conservative, int basis)
WeakOperators::ConvectionForm form, int basis)
{
size_t cmp = EV.size() / fe.basis(basis).size();
if (conservative) {
if (form == WeakOperators::CONSERVATIVE) {
for (size_t i = 1;i <= fe.basis(basis).size();i++)
for (size_t k = 1;k <= cmp;k++)
for (size_t l = 1;l <= cmp;l++)

View File

@ -18,6 +18,19 @@ class Vec3;
#include "FiniteElement.h"
#include "MatVec.h"
namespace WeakOperators
{
//! \brief Enum for the form of the convection term
enum ConvectionForm
{
CONVECTIVE = 0, //!< u_i u_j,i v_j
CONSERVATIVE = 1, //!< -u_i u_j v_i,j
SKEWSYMMETRIC = 2 //!< (u_i u_j,i v_j - u_i u_j v_i,j)/2
};
}
/*! \brief Common discrete operators using equal-ordered discretizations.
*/
@ -40,11 +53,12 @@ public:
//! \param[out] EM The element matrix to add contribution to
//! \param[in] fe The finite element to evaluate for
//! \param[in] U Advecting field
//! \param[in] conservative True to use the conservative formulation
//! \param[in] form Which form of the convective term to use
//! \param[in] basis Basis to use
static void Convection(Matrix& EM, const FiniteElement& fe,
const Vec3& U, const Tensor& dUdX, double scale,
bool conservative=false, int basis=1);
WeakOperators::ConvectionForm form=WeakOperators::CONVECTIVE,
int basis=1);
//! \brief Compute a divergence term.
//! \param[out] EM The element matrix to add contribution to
@ -135,11 +149,13 @@ public:
//! \param[in] dUdX Advected field gradient
//! \param[in] UC Advecting field
//! \param[in] scale Scaling factor for contribution
//! \param[in] conservative True to use the conservative formulation
//! \param[in] form Which form of the convective term to use
//! \param[in] basis Basis to use
static void Convection(Vector& EV, const FiniteElement& fe,
const Vec3& U, const Tensor& dUdX, const Vec3& UC,
double scale, bool conservative=false, int basis=1);
double scale,
WeakOperators::ConvectionForm form=WeakOperators::CONVECTIVE,
int basis=1);
//! \brief Compute a divergence term in a residual vector.
//! \param EV The element vector to add contribution to