https://mooseframework.inl.gov
NavierStokesLHDGAssemblyHelper.h
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 #pragma once
11 
13 
22 {
23 public:
25 
29  const TransientInterface * const ti,
30  const FEProblemBase & fe_problem,
31  SystemBase & sys,
32  const MooseMesh & mesh,
33  const THREAD_ID tid);
34 
35 protected:
41  const MooseArray<Number> & v_sol,
42  const unsigned int qp,
43  const unsigned int vel_component);
44 
50  const MooseArray<Number> & v_sol,
51  const unsigned int qp,
52  const unsigned int vel_component,
53  const unsigned int vel_j_component,
54  const MooseArray<std::vector<Real>> & phi,
55  const unsigned int j);
56 
64  void scalarVolumeResidual(const MooseArray<Gradient> & vel_gradient,
65  const unsigned int vel_component,
66  const Moose::Functor<Real> & body_force,
67  const MooseArray<Real> & JxW,
68  const libMesh::QBase & qrule,
69  const Elem * const current_elem,
70  const MooseArray<Point> & q_point,
71  DenseVector<Number> & scalar_re);
72 
83  void scalarVolumeJacobian(const unsigned int vel_component,
84  const MooseArray<Real> & JxW,
85  const libMesh::QBase & qrule,
86  DenseMatrix<Number> & scalar_vector_jac,
87  DenseMatrix<Number> & scalar_u_vel_jac,
88  DenseMatrix<Number> & scalar_v_vel_jac,
89  DenseMatrix<Number> & scalar_p_jac);
90 
98  void pressureVolumeResidual(const Moose::Functor<Real> & pressure_mms_forcing_function,
99  const MooseArray<Real> & JxW,
100  const libMesh::QBase & qrule,
101  const Elem * const current_elem,
102  const MooseArray<Point> & q_point,
103  DenseVector<Number> & pressure_re,
104  DenseVector<Number> & global_lm_re);
105 
116  void pressureVolumeJacobian(const MooseArray<Real> & JxW,
117  const libMesh::QBase & qrule,
118  DenseMatrix<Number> & p_u_vel_jac,
119  DenseMatrix<Number> & p_v_vel_jac,
120  DenseMatrix<Number> & p_global_lm_jac,
121  DenseMatrix<Number> & global_lm_p_jac);
122 
123  //
124  // Methods which are leveraged both on internal sides in the kernel and by the outflow bc
125  //
126 
127  void pressureFaceResidual(const MooseArray<Real> & JxW_face,
128  const libMesh::QBase & qrule_face,
129  const MooseArray<Point> & normals,
130  DenseVector<Number> & pressure_re);
131 
132  void pressureFaceJacobian(const MooseArray<Real> & JxW_face,
133  const libMesh::QBase & qrule_face,
134  const MooseArray<Point> & normals,
135  DenseMatrix<Number> & p_lm_u_vel_jac,
136  DenseMatrix<Number> & p_lm_v_vel_jac);
137 
138  void scalarFaceResidual(const MooseArray<Gradient> & vector_sol,
139  const MooseArray<Number> & scalar_sol,
140  const MooseArray<Number> & lm_sol,
141  const unsigned int vel_component,
142  const MooseArray<Real> & JxW_face,
143  const libMesh::QBase & qrule_face,
144  const MooseArray<Point> & normals,
145  DenseVector<Number> & scalar_re);
146 
147  void scalarFaceJacobian(const unsigned int vel_component,
148  const MooseArray<Real> & JxW_face,
149  const libMesh::QBase & qrule_face,
150  const MooseArray<Point> & normals,
151  DenseMatrix<Number> & scalar_vector_jac,
152  DenseMatrix<Number> & scalar_scalar_jac,
153  DenseMatrix<Number> & scalar_lm_jac,
154  DenseMatrix<Number> & scalar_p_jac,
155  DenseMatrix<Number> & scalar_lm_u_vel_jac,
156  DenseMatrix<Number> & scalar_lm_v_vel_jac);
157 
158  void lmFaceResidual(const MooseArray<Gradient> & vector_sol,
159  const MooseArray<Number> & scalar_sol,
160  const MooseArray<Number> & lm_sol,
161  const unsigned int vel_component,
162  const MooseArray<Real> & JxW_face,
163  const libMesh::QBase & qrule_face,
164  const MooseArray<Point> & normals,
165  const Elem * const neigh,
166  DenseVector<Number> & lm_re);
167 
168  void lmFaceJacobian(const unsigned int vel_component,
169  const MooseArray<Real> & JxW_face,
170  const libMesh::QBase & qrule_face,
171  const MooseArray<Point> & normals,
172  const Elem * const neigh,
173  DenseMatrix<Number> & lm_vec_jac,
174  DenseMatrix<Number> & lm_scalar_jac,
175  DenseMatrix<Number> & lm_lm_jac,
176  DenseMatrix<Number> & lm_p_jac,
177  DenseMatrix<Number> & lm_lm_u_vel_jac,
178  DenseMatrix<Number> & lm_lm_v_vel_jac);
179 
180  void pressureDirichletResidual(const std::array<const Moose::Functor<Real> *, 3> & dirichlet_vel,
181  const MooseArray<Real> & JxW_face,
182  const libMesh::QBase & qrule_face,
183  const MooseArray<Point> & normals,
184  const Elem * const current_elem,
185  const unsigned int current_side,
186  const MooseArray<Point> & q_point_face,
187  DenseVector<Number> & pressure_re);
188 
189  void scalarDirichletResidual(const MooseArray<Gradient> & vector_sol,
190  const MooseArray<Number> & scalar_sol,
191  const unsigned int vel_component,
192  const std::array<const Moose::Functor<Real> *, 3> & dirichlet_vel,
193  const MooseArray<Real> & JxW_face,
194  const libMesh::QBase & qrule_face,
195  const MooseArray<Point> & normals,
196  const Elem * const current_elem,
197  const unsigned int current_side,
198  const MooseArray<Point> & q_point_face,
199  DenseVector<Number> & scalar_re);
200 
201  void scalarDirichletJacobian(const unsigned int vel_component,
202  const MooseArray<Real> & JxW_face,
203  const libMesh::QBase & qrule_face,
204  const MooseArray<Point> & normals,
205  DenseMatrix<Number> & scalar_vector_jac,
206  DenseMatrix<Number> & scalar_scalar_jac,
207  DenseMatrix<Number> & scalar_pressure_jac);
208 
217 
219  const std::vector<dof_id_type> & _qv_dof_indices;
220  const std::vector<dof_id_type> & _v_dof_indices;
221  const std::vector<dof_id_type> & _lm_v_dof_indices;
222  const std::vector<dof_id_type> * const _qw_dof_indices;
223  const std::vector<dof_id_type> * const _w_dof_indices;
224  const std::vector<dof_id_type> * const _lm_w_dof_indices;
225  const std::vector<dof_id_type> & _p_dof_indices;
226  const std::vector<dof_id_type> * const _global_lm_dof_indices;
227 
233  const MooseArray<Number> * const _w_sol;
237 
239  const Real _rho;
240 
249 };
void pressureVolumeJacobian(const MooseArray< Real > &JxW, const libMesh::QBase &qrule, DenseMatrix< Number > &p_u_vel_jac, DenseMatrix< Number > &p_v_vel_jac, DenseMatrix< Number > &p_global_lm_jac, DenseMatrix< Number > &global_lm_p_jac)
Compute the volumetric contributions to the pressure Jacobian, e.g.
const std::vector< dof_id_type > *const _lm_w_dof_indices
const std::vector< dof_id_type > *const _global_lm_dof_indices
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
RealVectorValue rhoVelCrossVelJacobian(const MooseArray< Number > &u_sol, const MooseArray< Number > &v_sol, const unsigned int qp, const unsigned int vel_component, const unsigned int vel_j_component, const MooseArray< std::vector< Real >> &phi, const unsigned int j)
void pressureVolumeResidual(const Moose::Functor< Real > &pressure_mms_forcing_function, const MooseArray< Real > &JxW, const libMesh::QBase &qrule, const Elem *const current_elem, const MooseArray< Point > &q_point, DenseVector< Number > &pressure_re, DenseVector< Number > &global_lm_re)
Compute the volumetric contributions to the pressure residual, e.g.
const MooseVariableFE< Real > *const _w_var
const MooseVariableFE< Real > & _v_face_var
MeshBase & mesh
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)
const std::vector< dof_id_type > *const _qw_dof_indices
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.
void scalarVolumeJacobian(const unsigned int vel_component, const MooseArray< Real > &JxW, const libMesh::QBase &qrule, DenseMatrix< Number > &scalar_vector_jac, DenseMatrix< Number > &scalar_u_vel_jac, DenseMatrix< Number > &scalar_v_vel_jac, DenseMatrix< Number > &scalar_p_jac)
Compute the volumetric contributions to a velocity Jacobian.
RealVectorValue rhoVelCrossVelResidual(const MooseArray< Number > &u_sol, const MooseArray< Number > &v_sol, const unsigned int qp, const unsigned int vel_component)
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)
void scalarDirichletResidual(const MooseArray< Gradient > &vector_sol, const MooseArray< Number > &scalar_sol, const unsigned int vel_component, const std::array< const Moose::Functor< Real > *, 3 > &dirichlet_vel, const MooseArray< Real > &JxW_face, const libMesh::QBase &qrule_face, const MooseArray< Point > &normals, const Elem *const current_elem, const unsigned int current_side, const MooseArray< Point > &q_point_face, DenseVector< Number > &scalar_re)
const MooseVariableFE< RealVectorValue > & _grad_v_var
const std::vector< dof_id_type > & _lm_v_dof_indices
const MooseArray< Gradient > *const _qw_sol
const std::vector< dof_id_type > & _v_dof_indices
const MooseVariableFE< RealVectorValue > *const _grad_w_var
const std::vector< dof_id_type > *const _w_dof_indices
const MooseVariableFE< Real > & _pressure_var
const std::vector< dof_id_type > & _p_dof_indices
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void pressureDirichletResidual(const std::array< const Moose::Functor< Real > *, 3 > &dirichlet_vel, const MooseArray< Real > &JxW_face, const libMesh::QBase &qrule_face, const MooseArray< Point > &normals, const Elem *const current_elem, const unsigned int current_side, const MooseArray< Point > &q_point_face, DenseVector< Number > &pressure_re)
void pressureFaceResidual(const MooseArray< Real > &JxW_face, const libMesh::QBase &qrule_face, const MooseArray< Point > &normals, DenseVector< Number > &pressure_re)
NavierStokesLHDGAssemblyHelper(const MooseObject *moose_obj, MaterialPropertyInterface *mpi, MooseVariableDependencyInterface *mvdi, const TransientInterface *const ti, const FEProblemBase &fe_problem, SystemBase &sys, const MooseMesh &mesh, const THREAD_ID tid)
const MooseArray< Gradient > & _qv_sol
local solutions at quadrature points
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
const MooseArray< Number > *const _lm_w_sol
const MooseVariableFE< Real > *const _w_face_var
const MooseArray< Number > *const _w_sol
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)
void scalarDirichletJacobian(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_pressure_jac)
Implements all the methods for assembling a hybridized local discontinuous Galerkin (LDG-H)...
void scalarVolumeResidual(const MooseArray< Gradient > &vel_gradient, const unsigned int vel_component, const Moose::Functor< Real > &body_force, const MooseArray< Real > &JxW, const libMesh::QBase &qrule, const Elem *const current_elem, const MooseArray< Point > &q_point, DenseVector< Number > &scalar_re)
Compute the volumetric contributions to a velocity residual for a provided velocity gradient and stre...
unsigned int THREAD_ID
const MooseVariableScalar *const _enclosure_lm_var
const MooseArray< Number > *const _global_lm_dof_value