https://mooseframework.inl.gov
NavierStokesLHDGOutflowBC.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 
11 
13 
16 {
17  auto params = IntegratedBC::validParams();
19  params.addClassDescription("Implements an outflow boundary condition for use with a hybridized "
20  "discretization of the incompressible Navier-Stokes equations");
21  params.renameParam("variable", "u", "The x-component of velocity");
22  return params;
23 }
24 
26  : IntegratedBC(parameters),
27  NavierStokesLHDGAssemblyHelper(this, this, this, this, _fe_problem, _sys, _mesh, _tid),
28  _cached_side(libMesh::invalid_uint)
29 {
30 }
31 
32 void
34 {
35  checkCoupling();
36 }
37 
38 void
40 {
41  _cached_elem = nullptr;
43 }
44 
45 void
47 {
49  {
53  }
54 }
55 
56 void
58 {
59  const Elem * const neigh = _current_elem->neighbor_ptr(_current_side);
60 
67  _p_re.resize(_p_dof_indices.size());
68 
69  // qu, u, lm_u
73 
74  // qv, v, lm_v
78 
79  // p
81 
89 }
90 
91 void
93 {
94  const Elem * const neigh = _current_elem->neighbor_ptr(_current_side);
95 
120 
121  // qu, u, lm_u
124  _JxW,
125  *_qrule,
126  _normals,
128  _u_u_jac,
129  _u_lm_u_jac,
130  _u_p_jac,
131  _u_lm_u_jac,
132  _u_lm_v_jac);
133  lmFaceJacobian(0,
134  _JxW,
135  *_qrule,
136  _normals,
137  neigh,
139  _lm_u_u_jac,
141  _lm_u_p_jac,
144 
145  // qv, v, lm_v
148  _JxW,
149  *_qrule,
150  _normals,
152  _v_v_jac,
153  _v_lm_v_jac,
154  _v_p_jac,
155  _v_lm_u_jac,
156  _v_lm_v_jac);
157  lmFaceJacobian(1,
158  _JxW,
159  *_qrule,
160  _normals,
161  neigh,
163  _lm_v_v_jac,
165  _lm_v_p_jac,
168 
169  // p
171 
172  addJacobian(
179  addJacobian(
181  addJacobian(
183  addJacobian(
185  addJacobian(
187  addJacobian(
189  addJacobian(
196  addJacobian(
198  addJacobian(
200  addJacobian(
202  addJacobian(
204  addJacobian(
206  addJacobian(
208  addJacobian(
210 }
const std::vector< dof_id_type > & _lm_u_dof_indices
const MooseArray< Number > & _u_sol
void lmFaceJacobian(const unsigned int vel_component, const MooseArray< Real > &JxW_face, const libMesh::QBase &qrule_face, const MooseArray< Point > &normals, const Elem *const neigh, DenseMatrix< Number > &lm_vec_jac, DenseMatrix< Number > &lm_scalar_jac, DenseMatrix< Number > &lm_lm_jac, DenseMatrix< Number > &lm_p_jac, DenseMatrix< Number > &lm_lm_u_vel_jac, DenseMatrix< Number > &lm_lm_v_vel_jac)
const MooseVariableFE< Real > & _v_var
const unsigned int invalid_uint
const MooseArray< Point > & _normals
const std::vector< dof_id_type > & _u_dof_indices
void vectorFaceResidual(const MooseArray< Number > &lm_sol, const MooseArray< Real > &JxW_face, const libMesh::QBase &qrule_face, const MooseArray< Point > &normals, DenseVector< Number > &vector_re)
const Elem *const & _current_elem
void resize(const unsigned int n)
void addResiduals(Assembly &assembly, const Residuals &residuals, const Indices &dof_indices, Real scaling_factor)
virtual void computeOffDiagJacobian(unsigned int jvar) override
const MooseVariableFE< Real > & _v_face_var
static InputParameters validParams()
const MooseArray< libMesh::Gradient > & _qu_sol
void pressureFaceJacobian(const MooseArray< Real > &JxW_face, const libMesh::QBase &qrule_face, const MooseArray< Point > &normals, DenseMatrix< Number > &p_lm_u_vel_jac, DenseMatrix< Number > &p_lm_v_vel_jac)
void scalarFaceResidual(const MooseArray< Gradient > &vector_sol, const MooseArray< Number > &scalar_sol, const MooseArray< Number > &lm_sol, const unsigned int vel_component, const MooseArray< Real > &JxW_face, const libMesh::QBase &qrule_face, const MooseArray< Point > &normals, DenseVector< Number > &scalar_re)
const std::vector< dof_id_type > & _qv_dof_indices
Containers for dof indices.
The following methods are specializations for using the Parallel::packed_range_* routines for a vecto...
const std::vector< dof_id_type > & _qu_dof_indices
static InputParameters validParams()
void addJacobian(Assembly &assembly, const Residuals &residuals, const Indices &dof_indices, Real scaling_factor)
void scalarFaceJacobian(const unsigned int vel_component, const MooseArray< Real > &JxW_face, const libMesh::QBase &qrule_face, const MooseArray< Point > &normals, DenseMatrix< Number > &scalar_vector_jac, DenseMatrix< Number > &scalar_scalar_jac, DenseMatrix< Number > &scalar_lm_jac, DenseMatrix< Number > &scalar_p_jac, DenseMatrix< Number > &scalar_lm_u_vel_jac, DenseMatrix< Number > &scalar_lm_v_vel_jac)
NavierStokesLHDGOutflowBC(const InputParameters &parameters)
const MooseVariableFE< RealVectorValue > & _grad_v_var
const std::vector< dof_id_type > & _lm_v_dof_indices
const std::vector< dof_id_type > & _v_dof_indices
virtual void jacobianSetup() override
Assembly & _assembly
void vectorFaceJacobian(const MooseArray< Real > &JxW_face, const libMesh::QBase &qrule_face, const MooseArray< Point > &normals, DenseMatrix< Number > &vector_lm_jac)
virtual void computeResidual() override
const MooseVariableFE< Real > & _pressure_var
const MooseVariableFE< Real > & _u_face_var
const QBase *const & _qrule
const unsigned int & _current_side
const std::vector< dof_id_type > & _p_dof_indices
unsigned int _cached_side
A cache variable to prevent multiple computations of Jacobians.
void pressureFaceResidual(const MooseArray< Real > &JxW_face, const libMesh::QBase &qrule_face, const MooseArray< Point > &normals, DenseVector< Number > &pressure_re)
const MooseVariableFE< Real > & _u_var
const MooseArray< Number > & _lm_u_sol
void resize(const unsigned int new_m, const unsigned int new_n)
const MooseArray< Gradient > & _qv_sol
local solutions at quadrature points
virtual void initialSetup() override
const MooseVariableFE< RealVectorValue > & _grad_u_var
Implements an outflow boundary condition for use with a hybridized discretization of the incompressib...
void lmFaceResidual(const MooseArray< Gradient > &vector_sol, const MooseArray< Number > &scalar_sol, const MooseArray< Number > &lm_sol, const unsigned int vel_component, const MooseArray< Real > &JxW_face, const libMesh::QBase &qrule_face, const MooseArray< Point > &normals, const Elem *const neigh, DenseVector< Number > &lm_re)
Implements all the methods for assembling a hybridized local discontinuous Galerkin (LDG-H)...
virtual void computeJacobian() override
const MooseArray< Real > & _JxW
void scalingFactor(const std::vector< Real > &factor)