https://mooseframework.inl.gov
MassContinuityAssemblyHelper.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 "IPHDGAssemblyHelper.h"
13 #include "ADFunctorInterface.h"
14 
15 class MooseMesh;
16 
24 {
25 public:
27 
28  MassContinuityAssemblyHelper(const MooseObject * const moose_obj,
30  const TransientInterface * const ti,
31  const MooseMesh & mesh,
32  SystemBase & sys,
33  const Assembly & assembly,
34  const THREAD_ID tid,
35  const std::set<SubdomainID> & block_ids,
36  const std::set<BoundaryID> & boundary_ids);
37 
38  virtual void scalarVolume() override;
39  virtual void scalarFace() override;
40  virtual void lmFace() override;
41  virtual void scalarDirichlet(const Moose::Functor<Real> &) override
42  {
43  mooseError(
44  "This kernel is applied for the pressure variable, but the pressure does not itself appear "
45  "in the weak form this kernel implements, so the scalarDirichlet method does not fit");
46  }
47 
50 
52  const unsigned int _rz_radial_coord;
53 
55  std::vector<const ADVariableValue *> _interior_vels;
56 
58  std::vector<const ADVariableGradient *> _interior_vel_grads;
59 
61  std::vector<const Moose::Functor<ADReal> *> _face_vels;
62 };
void mooseError(Args &&... args)
const unsigned int _rz_radial_coord
The radial coordinate index for RZ coordinate systems.
Implements all the methods for assembling a hybridized interior penalty discontinuous Galerkin (IPDG-...
std::vector< const Moose::Functor< ADReal > * > _face_vels
The velocity functors used to evalute the velocity on element face quadrature points.
MeshBase & mesh
MassContinuityAssemblyHelper(const MooseObject *const moose_obj, MooseVariableDependencyInterface *const mvdi, const TransientInterface *const ti, const MooseMesh &mesh, SystemBase &sys, const Assembly &assembly, const THREAD_ID tid, const std::set< SubdomainID > &block_ids, const std::set< BoundaryID > &boundary_ids)
std::vector< const ADVariableGradient * > _interior_vel_grads
The velocity gradients at interior element quadrature points.
CoordinateSystemType
const Moose::CoordinateSystemType _coord_sys
The coordinate system.
virtual void scalarDirichlet(const Moose::Functor< Real > &) override
std::vector< const ADVariableValue * > _interior_vels
The velocities at interior element quadrature points.
unsigned int THREAD_ID