https://mooseframework.inl.gov
DiffusionLHDGAssemblyHelper.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 
12 #include "Moose.h"
13 #include "MooseTypes.h"
14 #include "MooseArray.h"
16 #include "ADFunctorInterface.h"
17 
18 #include "libmesh/vector_value.h"
19 #include <vector>
20 
21 template <typename>
22 class MooseVariableFE;
24 template <typename>
25 class MooseArray;
26 class SystemBase;
27 class MooseMesh;
28 class MooseObject;
31 class TransientInterface;
32 
41 {
42 public:
44 
45  DiffusionLHDGAssemblyHelper(const MooseObject * const moose_obj,
46  MaterialPropertyInterface * const mpi,
48  const TransientInterface * const ti,
49  const FEProblemBase & fe_problem,
50  SystemBase & sys,
51  const THREAD_ID tid);
52  void checkCoupling();
53 
54 protected:
61  void vectorVolumeResidual(const MooseArray<Gradient> & vector_sol,
62  const MooseArray<Number> & scalar_sol,
63  const MooseArray<Real> & JxW,
64  const libMesh::QBase & qrule,
65  DenseVector<Number> & vector_re);
66 
73  void vectorVolumeJacobian(const MooseArray<Real> & JxW,
74  const libMesh::QBase & qrule,
75  DenseMatrix<Number> & vector_vector_jac,
76  DenseMatrix<Number> & vector_scalar_jac);
77 
84  void scalarVolumeResidual(const MooseArray<Gradient> & vector_field,
85  const Moose::Functor<Real> & source,
86  const MooseArray<Real> & JxW,
87  const libMesh::QBase & qrule,
88  const Elem * const current_elem,
89  const MooseArray<Point> & q_point,
90  DenseVector<Number> & scalar_re);
91 
98  void scalarVolumeJacobian(const MooseArray<Real> & JxW,
99  const libMesh::QBase & qrule,
100  DenseMatrix<Number> & scalar_vector_jac);
101 
102  //
103  // Methods which can be leveraged both on internal sides in the kernel and by boundary conditions
104  //
105 
112  void vectorFaceResidual(const MooseArray<Number> & lm_sol,
113  const MooseArray<Real> & JxW_face,
114  const libMesh::QBase & qrule_face,
115  const MooseArray<Point> & normals,
116  DenseVector<Number> & vector_re);
117 
124  void vectorFaceJacobian(const MooseArray<Real> & JxW_face,
125  const libMesh::QBase & qrule_face,
126  const MooseArray<Point> & normals,
127  DenseMatrix<Number> & vector_lm_jac);
128 
133  void scalarFaceResidual(const MooseArray<Gradient> & vector_sol,
134  const MooseArray<Number> & scalar_sol,
135  const MooseArray<Number> & lm_sol,
136  const MooseArray<Real> & JxW_face,
137  const libMesh::QBase & qrule_face,
138  const MooseArray<Point> & normals,
139  DenseVector<Number> & scalar_re);
140 
145  void scalarFaceJacobian(const MooseArray<Real> & JxW_face,
146  const libMesh::QBase & qrule_face,
147  const MooseArray<Point> & normals,
148  DenseMatrix<Number> & scalar_vector_jac,
149  DenseMatrix<Number> & scalar_scalar_jac,
150  DenseMatrix<Number> & scalar_lm_jac);
151 
156  void lmFaceResidual(const MooseArray<Gradient> & vector_sol,
157  const MooseArray<Number> & scalar_sol,
158  const MooseArray<Number> & lm_sol,
159  const MooseArray<Real> & JxW_face,
160  const libMesh::QBase & qrule_face,
161  const MooseArray<Point> & normals,
162  DenseVector<Number> & lm_re);
163 
168  void lmFaceJacobian(const MooseArray<Real> & JxW_face,
169  const libMesh::QBase & qrule_face,
170  const MooseArray<Point> & normals,
171  DenseMatrix<Number> & lm_vec_jac,
172  DenseMatrix<Number> & lm_scalar_jac,
173  DenseMatrix<Number> & lm_lm_jac);
174 
178  void vectorDirichletResidual(const Moose::Functor<Real> & dirichlet_value,
179  const MooseArray<Real> & JxW_face,
180  const libMesh::QBase & qrule_face,
181  const MooseArray<Point> & normals,
182  const Elem * const current_elem,
183  const unsigned int current_side,
184  const MooseArray<Point> & q_point_face,
185  DenseVector<Number> & vector_re);
186 
190  void scalarDirichletResidual(const MooseArray<Gradient> & vector_sol,
191  const MooseArray<Number> & scalar_sol,
192  const Moose::Functor<Real> & dirichlet_value,
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 
205  void scalarDirichletJacobian(const MooseArray<Real> & JxW_face,
206  const libMesh::QBase & qrule_face,
207  const MooseArray<Point> & normals,
208  DenseMatrix<Number> & scalar_vector_jac,
209  DenseMatrix<Number> & scalar_scalar_jac);
210 
217  void createIdentityResidual(const MooseArray<Real> & JxW,
218  const libMesh::QBase & qrule,
219  const MooseArray<std::vector<Real>> & phi,
220  const MooseArray<Number> & sol,
221  DenseVector<Number> & re);
222 
226  void createIdentityJacobian(const MooseArray<Real> & JxW,
227  const libMesh::QBase & qrule,
228  const MooseArray<std::vector<Real>> & phi,
229  DenseMatrix<Number> & ke);
230 
234 
235  // Containers for dof indices
236  const std::vector<dof_id_type> & _qu_dof_indices;
237  const std::vector<dof_id_type> & _u_dof_indices;
238  const std::vector<dof_id_type> & _lm_u_dof_indices;
239 
240  // local solutions at quadrature points
244 
245  // Element shape functions
250 
251  // Face shape functions
255 
258 
261 
263  const Real _tau;
264 
266  const Elem * _cached_elem;
267 
268  // Local residual vectors
270 
271  // Local Jacobian matrices
275 
276 private:
279 
282 
285 };
const FEProblemBase & _dhah_fe_problem
A reference to the finite element problem used for coupling checks.
const std::vector< dof_id_type > & _lm_u_dof_indices
const MooseArray< std::vector< Real > > & _scalar_phi_face
const MooseArray< Number > & _u_sol
const MooseArray< std::vector< Real > > & _scalar_phi
void scalarDirichletResidual(const MooseArray< Gradient > &vector_sol, const MooseArray< Number > &scalar_sol, const Moose::Functor< Real > &dirichlet_value, 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)
Weakly imposes a Dirichlet condition for the scalar field in the scalar field equation.
const MooseArray< std::vector< Real > > & _div_vector_phi
void vectorVolumeJacobian(const MooseArray< Real > &JxW, const libMesh::QBase &qrule, DenseMatrix< Number > &vector_vector_jac, DenseMatrix< Number > &vector_scalar_jac)
Computes a local Jacobian matrix for the weak form: (q, v) + (u, div(v)) where q is the vector field ...
const std::vector< dof_id_type > & _u_dof_indices
Class for stuff related to variables.
Definition: Adaptivity.h:31
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)
Computes a local residual vector for the weak form: -<{u}, n*v> where {u} is the trace of the scalar ...
const Elem * _cached_elem
A data member used for determining when to compute the Jacobian.
const MooseArray< std::vector< Real > > & _lm_phi_face
This is a wrapper that forwards calls to the implementation, which can be switched out at any time wi...
const MooseArray< libMesh::Gradient > & _qu_sol
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
const MooseArray< std::vector< RealVectorValue > > & _grad_scalar_phi
Base class for a system (of equations)
Definition: SystemBase.h:84
Specialization of SubProblem for solving nonlinear equations plus auxiliary equations.
An interface for accessing Moose::Functors for systems that care about automatic differentiation, e.g.
const std::vector< dof_id_type > & _qu_dof_indices
Interface for objects that needs transient capabilities.
const SystemBase & _dhah_sys
A reference to the nonlinear system used for coupling checks.
void scalarDirichletJacobian(const MooseArray< Real > &JxW_face, const libMesh::QBase &qrule_face, const MooseArray< Point > &normals, DenseMatrix< Number > &scalar_vector_jac, DenseMatrix< Number > &scalar_scalar_jac)
Computes the Jacobian for a Dirichlet condition for the scalar field in the scalar field equation...
Every object that can be built by the factory should be derived from this class.
Definition: MooseObject.h:28
void vectorVolumeResidual(const MooseArray< Gradient > &vector_sol, const MooseArray< Number > &scalar_sol, const MooseArray< Real > &JxW, const libMesh::QBase &qrule, DenseVector< Number > &vector_re)
Computes a local residual vector for the weak form: (q, v) + (u, div(v)) where q is the vector field ...
static InputParameters validParams()
Implements all the methods for assembling a hybridized local discontinuous Galerkin (LDG-H)...
MooseMesh wraps a libMesh::Mesh object and enhances its capabilities by caching additional data and s...
Definition: MooseMesh.h:88
void scalarFaceJacobian(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)
Computes a local Jacobian matrix for the weak form: -<Dq*n, w> + < * (u - {u}) * n * n...
void lmFaceJacobian(const MooseArray< Real > &JxW_face, const libMesh::QBase &qrule_face, const MooseArray< Point > &normals, DenseMatrix< Number > &lm_vec_jac, DenseMatrix< Number > &lm_scalar_jac, DenseMatrix< Number > &lm_lm_jac)
Computes a local Jacobian matrix for the weak form: -<Dq*n, > + < * (u - {u}) * n * n...
void lmFaceResidual(const MooseArray< Gradient > &vector_sol, const MooseArray< Number > &scalar_sol, const MooseArray< Number > &lm_sol, const MooseArray< Real > &JxW_face, const libMesh::QBase &qrule_face, const MooseArray< Point > &normals, DenseVector< Number > &lm_re)
Computes a local residual vector for the weak form: -<Dq*n, > + < * (u - {u}) * n * n...
void vectorDirichletResidual(const Moose::Functor< Real > &dirichlet_value, 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 > &vector_re)
Weakly imposes a Dirichlet condition for the scalar field in the vector (gradient) equation...
const MaterialProperty< Real > & _diff
The diffusivity.
void scalarVolumeResidual(const MooseArray< Gradient > &vector_field, const Moose::Functor< Real > &source, const MooseArray< Real > &JxW, const libMesh::QBase &qrule, const Elem *const current_elem, const MooseArray< Point > &q_point, DenseVector< Number > &scalar_re)
Computes a local residual vector for the weak form: (Dq, grad(w)) - (f, w) where D is the diffusivity...
void vectorFaceJacobian(const MooseArray< Real > &JxW_face, const libMesh::QBase &qrule_face, const MooseArray< Point > &normals, DenseMatrix< Number > &vector_lm_jac)
Computes a local Jacobian matrix for the weak form: -<{u}, n*v> where {u} is the trace of the scalar ...
const TransientInterface & _ti
Reference to transient interface.
const MooseVariableFE< Real > & _u_face_var
forward declarations
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Real _tau
Our stabilization coefficient.
An interface for accessing Materials.
const MooseVariableFE< Real > & _u_var
const MooseObject & _moose_obj
A reference to our associated MooseObject for error reporting.
const MooseArray< Number > & _lm_u_sol
void createIdentityResidual(const MooseArray< Real > &JxW, const libMesh::QBase &qrule, const MooseArray< std::vector< Real >> &phi, const MooseArray< Number > &sol, DenseVector< Number > &re)
Creates residuals corresponding to the weak form (v, {u}), or stated simply this routine can be used ...
Class for scalar variables (they are different).
void scalarVolumeJacobian(const MooseArray< Real > &JxW, const libMesh::QBase &qrule, DenseMatrix< Number > &scalar_vector_jac)
Computes a local Jacobian matrix for the weak form: (Dq, grad(w)) - (f, w) where D is the diffusivity...
const MooseArray< std::vector< RealVectorValue > > & _vector_phi
const MooseVariableFE< RealVectorValue > & _grad_u_var
void scalarFaceResidual(const MooseArray< Gradient > &vector_sol, const MooseArray< Number > &scalar_sol, const MooseArray< Number > &lm_sol, const MooseArray< Real > &JxW_face, const libMesh::QBase &qrule_face, const MooseArray< Point > &normals, DenseVector< Number > &scalar_re)
Computes a local residual vector for the weak form: -<Dq*n, w> + < * (u - {u}) * n * n...
void createIdentityJacobian(const MooseArray< Real > &JxW, const libMesh::QBase &qrule, const MooseArray< std::vector< Real >> &phi, DenseMatrix< Number > &ke)
As above, but for the Jacobians.
const MooseArray< std::vector< RealVectorValue > > & _vector_phi_face
DiffusionLHDGAssemblyHelper(const MooseObject *const moose_obj, MaterialPropertyInterface *const mpi, MooseVariableDependencyInterface *const mvdi, const TransientInterface *const ti, const FEProblemBase &fe_problem, SystemBase &sys, const THREAD_ID tid)
unsigned int THREAD_ID
Definition: MooseTypes.h:209